CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types
PerigeeConversions Class Reference

#include <PerigeeConversions.h>

Public Member Functions

TrackCharge chargeFromPerigee (const PerigeeTrajectoryParameters &perigee) const
 
CurvilinearTrajectoryError curvilinearError (const PerigeeTrajectoryError &perigeeError, const GlobalTrajectoryParameters &gtp) const
 
PerigeeTrajectoryError ftsToPerigeeError (const FTS &originalFTS) const
 
PerigeeTrajectoryParameters ftsToPerigeeParameters (const FTS &originalFTS, const GlobalPoint &referencePoint, double &pt) const
 
AlgebraicMatrix55 jacobianCurvilinear2Perigee (const FreeTrajectoryState &fts) const
 
AlgebraicMatrix jacobianCurvilinear2Perigee_old (const FreeTrajectoryState &fts) const
 
AlgebraicMatrix66 jacobianParameters2Cartesian (const AlgebraicVector3 &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const
 
AlgebraicMatrix jacobianParameters2Cartesian_old (const AlgebraicVector &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const
 
AlgebraicMatrix55 jacobianPerigee2Curvilinear (const GlobalTrajectoryParameters &gtp) const
 
AlgebraicMatrix jacobianPerigee2Curvilinear_old (const GlobalTrajectoryParameters &gtp) const
 
GlobalVector momentumFromPerigee (const AlgebraicVector &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const
 
GlobalVector momentumFromPerigee (const AlgebraicVector3 &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const
 
GlobalVector momentumFromPerigee (const PerigeeTrajectoryParameters &parameters, double pt, const GlobalPoint &referencePoint) const
 
GlobalPoint positionFromPerigee (const PerigeeTrajectoryParameters &parameters, const GlobalPoint &referencePoint) const
 
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint (const AlgebraicVector &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicMatrix &theCovarianceMatrix, const MagneticField *field) const
 
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint (const AlgebraicVector3 &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicSymMatrix66 &theCovarianceMatrix, const MagneticField *field) const
 

Private Types

typedef FreeTrajectoryState FTS
 

Detailed Description

Class provides several methods to transform perigee parameters to and from various other parametrisations.

Definition at line 16 of file PerigeeConversions.h.

Member Typedef Documentation

Definition at line 18 of file PerigeeConversions.h.

Member Function Documentation

TrackCharge PerigeeConversions::chargeFromPerigee ( const PerigeeTrajectoryParameters perigee) const

This method returns the charge.

Definition at line 153 of file PerigeeConversions.cc.

References PerigeeTrajectoryParameters::charge().

154 {
155  return parameters.charge();
156 }
CurvilinearTrajectoryError PerigeeConversions::curvilinearError ( const PerigeeTrajectoryError perigeeError,
const GlobalTrajectoryParameters gtp 
) const

Definition at line 102 of file PerigeeConversions.cc.

References PerigeeTrajectoryError::covarianceMatrix().

Referenced by TrajectoryStateClosestToPoint::calculateFTS().

103 {
105  return CurvilinearTrajectoryError(ROOT::Math::Similarity(perigee2curv, perigeeError.covarianceMatrix()));
106 }
AlgebraicMatrix55 jacobianPerigee2Curvilinear(const GlobalTrajectoryParameters &gtp) const
const AlgebraicSymMatrix55 & covarianceMatrix() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
PerigeeTrajectoryError PerigeeConversions::ftsToPerigeeError ( const FTS originalFTS) const

Definition at line 63 of file PerigeeConversions.cc.

References FreeTrajectoryState::curvilinearError(), and CurvilinearTrajectoryError::matrix().

Referenced by MuonTrackingRegionBuilder::region(), and TrajectoryStateClosestToPoint::TrajectoryStateClosestToPoint().

64 {
65  AlgebraicSymMatrix55 errorMatrix = originalFTS.curvilinearError().matrix();
66  AlgebraicMatrix55 curv2perigee = jacobianCurvilinear2Perigee(originalFTS);
67  return PerigeeTrajectoryError(ROOT::Math::Similarity(curv2perigee,errorMatrix));
68 }
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const CurvilinearTrajectoryError & curvilinearError() const
AlgebraicMatrix55 jacobianCurvilinear2Perigee(const FreeTrajectoryState &fts) const
const AlgebraicSymMatrix55 & matrix() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
PerigeeTrajectoryParameters PerigeeConversions::ftsToPerigeeParameters ( const FTS originalFTS,
const GlobalPoint referencePoint,
double &  pt 
) const

This method calculates the perigee parameters from a given FTS and a reference point.

Definition at line 7 of file PerigeeConversions.cc.

References FreeTrajectoryState::charge(), epsilon, edm::hlt::Exception, MagneticField::inInverseGeV(), M_PI, GlobalTrajectoryParameters::magneticField(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), phi, PV3DBase< T, PVType, FrameType >::phi(), FreeTrajectoryState::position(), ExpressReco_HICollisions_FallBack::pt, mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::theta(), theta(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrajectoryStateClosestToPoint::TrajectoryStateClosestToPoint().

9 {
10  GlobalVector impactDistance = originalFTS.position() - referencePoint;
11 
12  pt = originalFTS.momentum().perp();
13  if (pt==0.) throw cms::Exception("PerigeeConversions", "Track with pt=0");
14 
15  double theta = originalFTS.momentum().theta();
16  double phi = originalFTS.momentum().phi();
17  double field = originalFTS.parameters().magneticField().inInverseGeV(originalFTS.position()).z();
18 // if (field==0.) throw cms::Exception("PerigeeConversions", "Field is 0") << " at " << originalFTS.position() << "\n" ;
19 
20  double positiveMomentumPhi = ( (phi>0) ? phi : (2*M_PI+phi) );
21  double positionPhi = impactDistance.phi();
22  double positivePositionPhi =
23  ( (positionPhi>0) ? positionPhi : (2*M_PI+positionPhi) );
24  double phiDiff = positiveMomentumPhi - positivePositionPhi;
25  if (phiDiff<0.0) phiDiff+= (2*M_PI);
26  double signEpsilon = ( (phiDiff > M_PI) ? -1.0 : 1.0);
27 
28  double epsilon = signEpsilon *
29  sqrt ( impactDistance.x()*impactDistance.x() +
30  impactDistance.y()*impactDistance.y() );
31 
32  // The track parameters:
33  AlgebraicVector5 theTrackParameters;
34 
35  double signTC = - originalFTS.charge();
36  bool isCharged = (signTC!=0) && (fabs(field)>1.e-10);
37  if (isCharged) {
38  theTrackParameters[0] = field / pt*signTC;
39  } else {
40  theTrackParameters[0] = 1 / pt;
41  }
42  theTrackParameters[1] = theta;
43  theTrackParameters[2] = phi;
44  theTrackParameters[3] = epsilon;
45  theTrackParameters[4] = impactDistance.z();
46  return PerigeeTrajectoryParameters(theTrackParameters, isCharged);
47 }
T perp() const
Definition: PV3DBase.h:66
const GlobalTrajectoryParameters & parameters() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
T y() const
Definition: PV3DBase.h:57
Geom::Theta< T > theta() const
TrackCharge charge() const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: DDAxes.h:10
T sqrt(T t)
Definition: SSEVec.h:28
T z() const
Definition: PV3DBase.h:58
GlobalVector momentum() const
GlobalPoint position() const
#define M_PI
Definition: BFit3D.cc:3
ROOT::Math::SVector< double, 5 > AlgebraicVector5
const MagneticField & magneticField() const
const double epsilon
T x() const
Definition: PV3DBase.h:56
Definition: DDAxes.h:10
AlgebraicMatrix55 PerigeeConversions::jacobianCurvilinear2Perigee ( const FreeTrajectoryState fts) const

Definition at line 226 of file PerigeeConversions.cc.

References alpha, funct::cos(), Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), ExpressReco_HICollisions_FallBack::e, Exhume::I, MagneticField::inInverseGeV(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), FreeTrajectoryState::momentum(), MultiGaussianStateTransform::N, L1TEmulatorMonitor_cff::p, FreeTrajectoryState::parameters(), FreeTrajectoryState::position(), FreeTrajectoryState::signedInverseMomentum(), funct::tan(), PV3DBase< T, PVType, FrameType >::theta(), FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), ExpressReco_HICollisions_FallBack::x, PV3DBase< T, PVType, FrameType >::y(), Gflash::Z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by jacobianCurvilinear2Perigee_old(), and PerigeeKinematicState::PerigeeKinematicState().

226  {
227 
228  GlobalVector p = fts.momentum();
229 
230  GlobalVector Z = GlobalVector(0.,0.,1.);
231  GlobalVector T = p.unit();
232  GlobalVector U = Z.cross(T).unit();;
233  GlobalVector V = T.cross(U);
234 
235  GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.); //opposite to track dir.
236  I = I.unit();
237  GlobalVector J(-I.y(), I.x(),0.); //counterclockwise rotation
238  GlobalVector K(Z);
239  GlobalPoint x = fts.position();
241  GlobalVector H = B.unit();
242  GlobalVector HxT = H.cross(T);
243  GlobalVector N = HxT.unit();
244  double alpha = HxT.mag();
245  double qbp = fts.signedInverseMomentum();
246  double Q = -B.mag() * qbp;
247  double alphaQ = alpha * Q;
248 
249  double lambda = 0.5 * M_PI - p.theta();
250  double coslambda = cos(lambda), tanlambda = tan(lambda);
251 
252  double TI = T.dot(I);
253  double NU = N.dot(U);
254  double NV = N.dot(V);
255  double UI = U.dot(I);
256  double VI = V.dot(I);
257  double UJ = U.dot(J);
258  double VJ = V.dot(J);
259  double UK = U.dot(K);
260  double VK = V.dot(K);
261 
262  AlgebraicMatrix55 jac;
263 
264  if( fabs(fts.transverseCurvature())<1.e-10 ) {
265  jac(0,0) = 1/coslambda;
266  jac(0,1) = tanlambda/coslambda/fts.parameters().momentum().mag();
267  }else{
268  double Bz = B.z();
269  jac(0,0) = -Bz/coslambda;
270  jac(0,1) = -Bz * tanlambda/coslambda*qbp;
271  jac(1,3) = alphaQ * NV * UI/TI;
272  jac(1,4) = alphaQ * NV * VI/TI;
273  jac(0,3) = -jac(0,1) * jac(1,3);
274  jac(0,4) = -jac(0,1) * jac(1,4);
275  jac(2,3) = -alphaQ/coslambda * NU * UI/TI;
276  jac(2,4) = -alphaQ/coslambda * NU * VI/TI;
277  }
278  jac(1,1) = -1.;
279  jac(2,2) = 1.;
280  jac(3,3) = VK/TI;
281  jac(3,4) = -UK/TI;
282  jac(4,3) = -VJ/TI;
283  jac(4,4) = UJ/TI;
284 
285  return jac;
286 
287 }
const double Z[kNumberCalorimeter]
float alpha
Definition: AMPTWrapper.h:95
const GlobalTrajectoryParameters & parameters() const
T y() const
Definition: PV3DBase.h:57
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
T mag() const
Definition: PV3DBase.h:61
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
const std::complex< double > I
Definition: I.h:8
GlobalVector momentum() const
Vector3DBase unit() const
Definition: Vector3DBase.h:57
GlobalPoint position() const
#define M_PI
Definition: BFit3D.cc:3
double transverseCurvature() const
const MagneticField & magneticField() const
T x() const
Definition: PV3DBase.h:56
double signedInverseMomentum() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
Global3DVector GlobalVector
Definition: GlobalVector.h:10
AlgebraicMatrix PerigeeConversions::jacobianCurvilinear2Perigee_old ( const FreeTrajectoryState fts) const

Jacobians of tranformations between curvilinear frame at point of closest approach in transverse plane and perigee frame. The fts must therefore be given at exactly this point in order to yield the correct Jacobians.

Definition at line 220 of file PerigeeConversions.cc.

References asHepMatrix(), and jacobianCurvilinear2Perigee().

220  {
222 }
CLHEP::HepMatrix asHepMatrix(const ROOT::Math::SMatrix< double, N1, N2, typename ROOT::Math::MatRepStd< double, N1, N2 > > &rm)
Definition: Migration.h:49
AlgebraicMatrix55 jacobianCurvilinear2Perigee(const FreeTrajectoryState &fts) const
AlgebraicMatrix66 PerigeeConversions::jacobianParameters2Cartesian ( const AlgebraicVector3 momentum,
const GlobalPoint position,
const TrackCharge charge,
const MagneticField field 
) const

Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian

Definition at line 193 of file PerigeeConversions.cc.

References ecalTB2006H4_GenSimDigiReco_cfg::bField, DeDxDiscriminatorTools::charge(), funct::cos(), ExpressReco_HICollisions_FallBack::e, edm::hlt::Exception, MagneticField::inInverseGeV(), funct::sin(), funct::tan(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by jacobianParameters2Cartesian_old(), and KinematicPerigeeConversions::jacobianParameters2Kinematic().

196 {
197  if (momentum[0]==0.) throw cms::Exception("PerigeeConversions", "Track with rho=0");
198  double factor = 1.;
199  double bField = field->inInverseGeV(position).z();
200  if (charge!=0 && fabs(bField)>1.e-10) {
201 // if (bField==0.) throw cms::Exception("PerigeeConversions", "Field is 0");
202  factor = -bField*charge;
203  }
204  AlgebraicMatrix66 frameTransJ;
205  frameTransJ(0,0) = 1;
206  frameTransJ(1,1) = 1;
207  frameTransJ(2,2) = 1;
208  frameTransJ(3,3) = - factor * cos(momentum[2]) / (momentum[0]*momentum[0]);
209  frameTransJ(4,3) = - factor * sin(momentum[2]) / (momentum[0]*momentum[0]);
210  frameTransJ(5,3) = - factor / tan(momentum[1]) / (momentum[0]*momentum[0]);
211 
212  frameTransJ(3,5) = - factor * sin(momentum[2]) / (momentum[0]);
213  frameTransJ(4,5) = factor * cos(momentum[2]) / (momentum[0]);
214  frameTransJ(5,4) = - factor / (momentum[0]*sin(momentum[1])*sin(momentum[1]));
215 
216  return frameTransJ;
217 }
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepStd< double, 6, 6 > > AlgebraicMatrix66
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double charge(const std::vector< uint8_t > &Ampls)
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
AlgebraicMatrix PerigeeConversions::jacobianParameters2Cartesian_old ( const AlgebraicVector momentum,
const GlobalPoint position,
const TrackCharge charge,
const MagneticField field 
) const

Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian

Definition at line 186 of file PerigeeConversions.cc.

References asHepMatrix(), and jacobianParameters2Cartesian().

188  {
189  return asHepMatrix(jacobianParameters2Cartesian(asSVector<3>(momentum), position, charge, field));
190 }
CLHEP::HepMatrix asHepMatrix(const ROOT::Math::SMatrix< double, N1, N2, typename ROOT::Math::MatRepStd< double, N1, N2 > > &rm)
Definition: Migration.h:49
double charge(const std::vector< uint8_t > &Ampls)
AlgebraicMatrix66 jacobianParameters2Cartesian(const AlgebraicVector3 &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const
AlgebraicMatrix55 PerigeeConversions::jacobianPerigee2Curvilinear ( const GlobalTrajectoryParameters gtp) const

Definition at line 296 of file PerigeeConversions.cc.

References alpha, funct::cos(), Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), ExpressReco_HICollisions_FallBack::e, Exhume::I, MagneticField::inInverseGeV(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), MultiGaussianStateTransform::N, L1TEmulatorMonitor_cff::p, GlobalTrajectoryParameters::position(), GlobalTrajectoryParameters::signedInverseMomentum(), funct::sin(), PV3DBase< T, PVType, FrameType >::theta(), GlobalTrajectoryParameters::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), ExpressReco_HICollisions_FallBack::x, PV3DBase< T, PVType, FrameType >::y(), Gflash::Z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by jacobianPerigee2Curvilinear_old().

296  {
297 
298  GlobalVector p = gtp.momentum();
299 
300  GlobalVector Z = GlobalVector(0.,0.,1.);
301  GlobalVector T = p.unit();
302  GlobalVector U = Z.cross(T).unit();;
303  GlobalVector V = T.cross(U);
304 
305  GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.); //opposite to track dir.
306  I = I.unit();
307  GlobalVector J(-I.y(), I.x(),0.); //counterclockwise rotation
308  GlobalVector K(Z);
309  GlobalPoint x = gtp.position();
311  GlobalVector H = B.unit();
312  GlobalVector HxT = H.cross(T);
313  GlobalVector N = HxT.unit();
314  double alpha = HxT.mag();
315  double qbp = gtp.signedInverseMomentum();
316  double Q = -B.mag() * qbp;
317  double alphaQ = alpha * Q;
318 
319  double lambda = 0.5 * M_PI - p.theta();
320  double coslambda = cos(lambda), sinlambda = sin(lambda);
321  double mqbpt = -1./coslambda * qbp;
322 
323  double TJ = T.dot(J);
324  double TK = T.dot(K);
325  double NU = N.dot(U);
326  double NV = N.dot(V);
327  double UJ = U.dot(J);
328  double VJ = V.dot(J);
329  double UK = U.dot(K);
330  double VK = V.dot(K);
331 
332  AlgebraicMatrix55 jac;
333 
334  if( fabs(gtp.transverseCurvature())<1.e-10 ) {
335  jac(0,0) = coslambda;
336  jac(0,1) = sinlambda/coslambda/gtp.momentum().mag();
337  }else{
338  jac(0,0) = -coslambda/B.z();
339  jac(0,1) = -sinlambda * mqbpt;
340  jac(1,3) = -alphaQ * NV * TJ;
341  jac(1,4) = -alphaQ * NV * TK;
342  jac(2,3) = -alphaQ/coslambda * NU * TJ;
343  jac(2,4) = -alphaQ/coslambda * NU * TK;
344  }
345  jac(1,1) = -1.;
346  jac(2,2) = 1.;
347  jac(3,3) = UJ;
348  jac(3,4) = UK;
349  jac(4,3) = VJ;
350  jac(4,4) = VK;
351 
352  return jac;
353 }
const double Z[kNumberCalorimeter]
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:57
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
T mag() const
Definition: PV3DBase.h:61
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const std::complex< double > I
Definition: I.h:8
Vector3DBase unit() const
Definition: Vector3DBase.h:57
#define M_PI
Definition: BFit3D.cc:3
const MagneticField & magneticField() const
T x() const
Definition: PV3DBase.h:56
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
Global3DVector GlobalVector
Definition: GlobalVector.h:10
AlgebraicMatrix PerigeeConversions::jacobianPerigee2Curvilinear_old ( const GlobalTrajectoryParameters gtp) const

Definition at line 291 of file PerigeeConversions.cc.

References asHepMatrix(), and jacobianPerigee2Curvilinear().

291  {
293 }
CLHEP::HepMatrix asHepMatrix(const ROOT::Math::SMatrix< double, N1, N2, typename ROOT::Math::MatRepStd< double, N1, N2 > > &rm)
Definition: Migration.h:49
AlgebraicMatrix55 jacobianPerigee2Curvilinear(const GlobalTrajectoryParameters &gtp) const
GlobalVector PerigeeConversions::momentumFromPerigee ( const AlgebraicVector momentum,
const TrackCharge charge,
const GlobalPoint referencePoint,
const MagneticField field 
) const

This method returns the (Cartesian) momentum. The parameters need not be the full perigee parameters, as long as the first 3 parameters are the transverse curvature, theta and phi.

Definition at line 127 of file PerigeeConversions.cc.

Referenced by TrajectoryStateClosestToPoint::calculateFTS(), and TrajectoryStateClosestToPoint::momentum().

128  {
129  return momentumFromPerigee(asSVector<3>(momentum), charge, referencePoint, field);
130  }
double charge(const std::vector< uint8_t > &Ampls)
GlobalVector momentumFromPerigee(const AlgebraicVector &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const
GlobalVector PerigeeConversions::momentumFromPerigee ( const AlgebraicVector3 momentum,
const TrackCharge charge,
const GlobalPoint referencePoint,
const MagneticField field 
) const

Definition at line 133 of file PerigeeConversions.cc.

References abs, funct::cos(), ExpressReco_HICollisions_FallBack::e, edm::hlt::Exception, MagneticField::inInverseGeV(), ExpressReco_HICollisions_FallBack::pt, funct::sin(), funct::tan(), and PV3DBase< T, PVType, FrameType >::z().

135 {
136  double pt;
137  if (momentum[0]==0.) throw cms::Exception("PerigeeConversions", "Track with rho=0");
138 
139  double bz = fabs(field->inInverseGeV(referencePoint).z());
140  if ( charge!=0 && bz>1.e-10 ) {
141  pt = std::abs(bz/momentum[0]);
142  if (pt<1.e-10) throw cms::Exception("PerigeeConversions", "pt is 0");
143  } else {
144  pt = 1 / momentum[0];
145  }
146 
147  return GlobalVector(cos(momentum[2]) * pt,
148  sin(momentum[2]) * pt,
149  pt/tan(momentum[1]));
150 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define abs(x)
Definition: mlp_lapack.h:159
double charge(const std::vector< uint8_t > &Ampls)
virtual GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalVector PerigeeConversions::momentumFromPerigee ( const PerigeeTrajectoryParameters parameters,
double  pt,
const GlobalPoint referencePoint 
) const

This method returns the (Cartesian) momentum from the PerigeeTrajectoryParameters

Definition at line 119 of file PerigeeConversions.cc.

References funct::cos(), PerigeeTrajectoryParameters::phi(), funct::sin(), funct::tan(), and PerigeeTrajectoryParameters::theta().

120 {
121  return GlobalVector(cos(parameters.phi()) * pt,
122  sin(parameters.phi()) * pt,
123  pt / tan(parameters.theta()));
124 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalPoint PerigeeConversions::positionFromPerigee ( const PerigeeTrajectoryParameters parameters,
const GlobalPoint referencePoint 
) const

This method returns the position (on the helix) at which the parameters are defined

Definition at line 109 of file PerigeeConversions.cc.

References funct::cos(), funct::sin(), PerigeeTrajectoryParameters::vector(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrajectoryStateClosestToPoint::calculateFTS(), and TrajectoryStateClosestToPoint::position().

110 {
111  AlgebraicVector5 theVector = parameters.vector();
112  return GlobalPoint(theVector[3]*sin(theVector[2])+referencePoint.x(),
113  -theVector[3]*cos(theVector[2])+referencePoint.y(),
114  theVector[4]+referencePoint.z());
115 }
const AlgebraicVector5 & vector() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
ROOT::Math::SVector< double, 5 > AlgebraicVector5
T x() const
Definition: PV3DBase.h:56
TrajectoryStateClosestToPoint PerigeeConversions::trajectoryStateClosestToPoint ( const AlgebraicVector momentum,
const GlobalPoint referencePoint,
const TrackCharge charge,
const AlgebraicMatrix theCovarianceMatrix,
const MagneticField field 
) const

Public constructor. This constructor takes a momentum, with parameters (transverse curvature, theta, phi) and a position, which is both the reference position and the position at which the momentum is defined. The covariance matrix is defined for these 6 parameters, in the order (x, y, z, transverse curvature, theta, phi).

Parameters
fieldFIXME !!! why not Sym !!??

Definition at line 159 of file PerigeeConversions.cc.

Referenced by PerigeeMultiLTS::createRefittedTrackState(), and PerigeeLinearizedTrackState::createRefittedTrackState().

161  {
162  AlgebraicSymMatrix sym; sym.assign(theCovarianceMatrix); // below, this was used for Matrix => SymMatrix
163  return trajectoryStateClosestToPoint(asSVector<3>(momentum), referencePoint,
164  charge, asSMatrix<6>(sym), field);
165 
166  }
double charge(const std::vector< uint8_t > &Ampls)
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const AlgebraicVector &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicMatrix &theCovarianceMatrix, const MagneticField *field) const
CLHEP::HepSymMatrix AlgebraicSymMatrix
TrajectoryStateClosestToPoint PerigeeConversions::trajectoryStateClosestToPoint ( const AlgebraicVector3 momentum,
const GlobalPoint referencePoint,
const TrackCharge charge,
const AlgebraicSymMatrix66 theCovarianceMatrix,
const MagneticField field 
) const

Definition at line 170 of file PerigeeConversions.cc.

173 {
175  (momentum, referencePoint, charge, field);
176  CartesianTrajectoryError cartesianTrajErr(ROOT::Math::Similarity(param2cart, theCovarianceMatrix));
177 
178  FTS theFTS(GlobalTrajectoryParameters(referencePoint,
179  momentumFromPerigee(momentum, charge, referencePoint, field), charge,
180  field), cartesianTrajErr);
181 
182  return TrajectoryStateClosestToPoint(theFTS, referencePoint);
183 }
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepStd< double, 6, 6 > > AlgebraicMatrix66
double charge(const std::vector< uint8_t > &Ampls)
AlgebraicMatrix66 jacobianParameters2Cartesian(const AlgebraicVector3 &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const
GlobalVector momentumFromPerigee(const AlgebraicVector &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const