CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TSCPBuilderNoMaterial.cc
Go to the documentation of this file.
11 
14  const GlobalPoint& referencePoint) const
15 {
16  if (positionEqual(referencePoint, originalFTS.position()))
17  return constructTSCP(originalFTS, referencePoint);
18 
19  // Now do the propagation or whatever...
20 
21  PairBoolFTS newStatePair =
22  createFTSatTransverseImpactPoint(originalFTS, referencePoint);
23  if (newStatePair.first) {
24  return constructTSCP(newStatePair.second, referencePoint);
25  } else {
27  }
28 }
29 
32  const GlobalPoint& referencePoint) const
33 {
34  if (positionEqual(referencePoint, originalTSOS.globalPosition()))
35  return constructTSCP(*originalTSOS.freeState(), referencePoint);
36 
37  // Now do the propagation
38 
39  PairBoolFTS newStatePair =
40  createFTSatTransverseImpactPoint(*originalTSOS.freeState(), referencePoint);
41  if (newStatePair.first) {
42  return constructTSCP(newStatePair.second, referencePoint);
43  } else {
45  }
46 }
47 
50  const FTS& originalFTS, const GlobalPoint& referencePoint) const
51 {
52  //
53  // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um
54  // difference in transversal position at 10m.
55  //
56  if( fabs(originalFTS.transverseCurvature())<1.e-10 ) {
57  return createFTSatTransverseImpactPointNeutral(originalFTS, referencePoint);
58  } else {
59  return createFTSatTransverseImpactPointCharged(originalFTS, referencePoint);
60  }
61 }
62 
65  const FTS& originalFTS, const GlobalPoint& referencePoint) const
66 {
67 
68  GlobalVector pvecOrig = originalFTS.momentum();
69  GlobalPoint xvecOrig = originalFTS.position();
70  double kappa = originalFTS.transverseCurvature();
71  double pxOrig = pvecOrig.x();
72  double pyOrig = pvecOrig.y();
73  double pzOrig = pvecOrig.z();
74  double xOrig = xvecOrig.x();
75  double yOrig = xvecOrig.y();
76  double zOrig = xvecOrig.z();
77 
78 // double fac = 1./originalFTS.charge()/MagneticField::inInverseGeV(referencePoint).z();
79  double fac = 1./originalFTS.charge()/
80  (originalFTS.parameters().magneticField().inInverseGeV(referencePoint).z());
81  GlobalVectorDouble xOrig2Centre = GlobalVectorDouble(fac * pyOrig, -fac * pxOrig, 0.);
82  GlobalVectorDouble xOrigProj = GlobalVectorDouble(xOrig, yOrig, 0.);
83  GlobalVectorDouble xRefProj = GlobalVectorDouble(referencePoint.x(), referencePoint.y(), 0.);
84  GlobalVectorDouble deltax = xRefProj-xOrigProj-xOrig2Centre;
85  GlobalVectorDouble ndeltax = deltax.unit();
86 
88  Surface::PositionType pos(referencePoint);
89  // Need to define plane with orientation as the
90  // ImpactPointSurface
91  GlobalVector X(ndeltax.x(), ndeltax.y(), ndeltax.z());
92  GlobalVector Y(0.,0.,1.);
93  Surface::RotationType rot(X,Y);
94  BoundPlane* plane = new BoundPlane(pos,rot);
95  // Using Teddy's HelixBarrelPlaneCrossingByCircle for general barrel planes.
96  // A large variety of other,
97  // direct solutions turned out to be not so stable numerically.
99  planeCrossing(HelixPlaneCrossing::PositionType(xOrig, yOrig, zOrig),
100  HelixPlaneCrossing::DirectionType(pxOrig, pyOrig, pzOrig),
101  kappa, direction);
102  std::pair<bool,double> propResult = planeCrossing.pathLength(*plane);
103  if ( !propResult.first ) {
104  edm::LogWarning ("TSCPBuilderNoMaterial") << "Propagation to perigee plane failed!";
105  return PairBoolFTS(false, FreeTrajectoryState() );
106  }
107  double s = propResult.second;
108  HelixPlaneCrossing::PositionType xGen = planeCrossing.position(s);
109  GlobalPoint xPerigee = GlobalPoint(xGen.x(),xGen.y(),xGen.z());
110  // direction (reconverted to GlobalVector, renormalised)
111  HelixPlaneCrossing::DirectionType pGen = planeCrossing.direction(s);
112  pGen *= pvecOrig.mag()/pGen.mag();
113  GlobalVector pPerigee = GlobalVector(pGen.x(),pGen.y(),pGen.z());
114  delete plane;
115 
116  if (originalFTS.hasError()) {
117  const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
118  AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
119  pPerigee, s);
120  const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
121  CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
122 
123  return PairBoolFTS(true,
124  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
125  &(originalFTS.parameters().magneticField())), cte) );
126  }
127  else {
128  return PairBoolFTS(true,
129  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
130  &(originalFTS.parameters().magneticField()))) );
131  }
132 
133 }
134 
135 
138  const GlobalPoint& referencePoint) const
139 {
140 
141  GlobalPoint xvecOrig = originalFTS.position();
142  double phi = originalFTS.momentum().phi();
143  double theta = originalFTS.momentum().theta();
144  double xOrig = xvecOrig.x();
145  double yOrig = xvecOrig.y();
146  double zOrig = xvecOrig.z();
147  double xR = referencePoint.x();
148  double yR = referencePoint.y();
149 
150  double s2D = (xR - xOrig) * cos(phi) + (yR - yOrig) * sin(phi);
151  double s = s2D / sin(theta);
152  double xGen = xOrig + s2D*cos(phi);
153  double yGen = yOrig + s2D*sin(phi);
154  double zGen = zOrig + s2D/tan(theta);
155  GlobalPoint xPerigee = GlobalPoint(xGen, yGen, zGen);
156 
157  GlobalVector pPerigee = originalFTS.momentum();
158 
159  if (originalFTS.hasError()) {
160  const AlgebraicSymMatrix55 &errorMatrix = originalFTS.curvilinearError().matrix();
161  AnalyticalCurvilinearJacobian curvilinJacobian(originalFTS.parameters(), xPerigee,
162  pPerigee, s);
163  const AlgebraicMatrix55 &jacobian = curvilinJacobian.jacobian();
164  CurvilinearTrajectoryError cte( ROOT::Math::Similarity(jacobian, errorMatrix) );
165 
166  return PairBoolFTS(true,
167  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
168  &(originalFTS.parameters().magneticField())), cte));
169  }
170  else {
171  return PairBoolFTS(true,
172  FreeTrajectoryState(GlobalTrajectoryParameters(xPerigee, pPerigee, originalFTS.charge(),
173  &(originalFTS.parameters().magneticField()))) );
174  }
175 
176 }
virtual PositionType position(double s) const
T y() const
Cartesian y coordinate.
T x() const
Cartesian x coordinate.
virtual TrajectoryStateClosestToPoint operator()(const FTS &originalFTS, const GlobalPoint &referencePoint) const
const GlobalTrajectoryParameters & parameters() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
TrajectoryStateClosestToPoint constructTSCP(const FTS &originalFTS, const GlobalPoint &referencePoint) const
PairBoolFTS createFTSatTransverseImpactPointCharged(const FTS &originalFTS, const GlobalPoint &referencePoint) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
#define X(str)
Definition: MuonsGrabber.cc:49
GlobalPoint globalPosition() const
Geom::Theta< T > theta() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
PropagationDirection
TrackCharge charge() const
std::pair< bool, FreeTrajectoryState > PairBoolFTS
const CurvilinearTrajectoryError & curvilinearError() 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.
T mag() const
Definition: PV3DBase.h:61
FreeTrajectoryState * freeState(bool withErrors=true) const
T z() const
Cartesian z coordinate.
bool positionEqual(const GlobalPoint &ptB, const GlobalPoint &ptA) const
T z() const
Definition: PV3DBase.h:58
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Vector3DBase< double, GlobalTag > GlobalVectorDouble
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
virtual const AlgebraicMatrix55 & jacobian() const
GlobalVector momentum() const
virtual std::pair< bool, double > pathLength(const Plane &)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
GlobalPoint position() const
double transverseCurvature() const
const AlgebraicSymMatrix55 & matrix() const
const MagneticField & magneticField() const
PairBoolFTS createFTSatTransverseImpactPointNeutral(const FTS &originalFTS, const GlobalPoint &referencePoint) const
string s
Definition: asciidump.py:422
T x() const
Definition: PV3DBase.h:56
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
virtual DirectionType direction(double s) const
Global3DVector GlobalVector
Definition: GlobalVector.h:10
PairBoolFTS createFTSatTransverseImpactPoint(const FTS &originalFTS, const GlobalPoint &referencePoint) const
Definition: DDAxes.h:10