CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PerigeeLinearizedTrackState.cc
Go to the documentation of this file.
7 
8 
9 
11 {
12  GlobalPoint paramPt(theLinPoint);
13 
14  thePredState = builder(theTSOS, paramPt);
15  if unlikely(!thePredState.isValid()) return;
16 
17  double field = theTrack.field()->inInverseGeV(thePredState.theState().position()).z();
18 
19  if ((std::abs(theCharge)<1e-5)||(fabs(field)<1.e-10)){
20  //neutral track
22  } else {
23  //charged track
25  }
26 
27  jacobiansAvailable = true;
28 }
29 
30 
31 
33 {
34  const PerigeeLinearizedTrackState* otherP =
35  dynamic_cast<const PerigeeLinearizedTrackState*>(&other);
36  if (otherP == 0) {
37  throw VertexException("PerigeeLinearizedTrackState: don't know how to compare myself to non-perigee track state");
38  }
39  return (otherP->track() == theTrack);
40 }
41 
42 
44 {
45  const PerigeeLinearizedTrackState* otherP =
46  dynamic_cast<const PerigeeLinearizedTrackState*>(other.get());
47  if (otherP == 0) {
48  throw VertexException("PerigeeLinearizedTrackState: don't know how to compare myself to non-perigee track state");
49  }
50  return (otherP->track() == theTrack);
51 }
52 
53 
56  (const GlobalPoint & newLP) const
57 {
59  new PerigeeLinearizedTrackState(newLP, track(), theTSOS));
60 }
61 
64  const GlobalPoint & vertexPosition,
65  const AlgebraicVector3 & vectorParameters,
66  const AlgebraicSymMatrix66 & covarianceMatrix) const
67 {
68  TrajectoryStateClosestToPoint refittedTSCP =
70  vectorParameters, vertexPosition, charge(), covarianceMatrix, theTrack.field());
71  return RefCountedRefittedTrackState(new PerigeeRefittedTrackState(refittedTSCP, vectorParameters));
72 }
73 
74 std::vector< PerigeeLinearizedTrackState::RefCountedLinearizedTrackState >
76 {
77  std::vector<RefCountedLinearizedTrackState> result; result.reserve(1);
78  result.push_back(RefCountedLinearizedTrackState(
79  const_cast<PerigeeLinearizedTrackState*>(this)));
80  return result;
81 }
82 
83 
85  const RefCountedRefittedTrackState & theRefittedState) const
86 {
87  auto p = theRefittedState->position();
88  AlgebraicVector3 vertexPosition(p.x(),p.y(),p.z());
89  AlgebraicVector3 momentum = theRefittedState->momentumVector();
90  if ((momentum(2)*predictedStateMomentumParameters()(2) < 0)&&(fabs(momentum(2))>M_PI/2) ) {
91  if (predictedStateMomentumParameters()(2) < 0.) momentum(2)-= 2*M_PI;
92  if (predictedStateMomentumParameters()(2) > 0.) momentum(2)+= 2*M_PI;
93  }
94  AlgebraicVectorN param = constantTerm() +
95  positionJacobian() * vertexPosition +
96  momentumJacobian() * momentum;
97  if (param(2) > M_PI) param(2)-= 2*M_PI;
98  if (param(2) < -M_PI) param(2)+= 2*M_PI;
99 
100  return param;
101 }
102 
103 
105 {
106  if (parameters(2) > M_PI) parameters(2)-= 2*M_PI;
107  if (parameters(2) < -M_PI) parameters(2)+= 2*M_PI;
108 }
109 
111 {
112  GlobalPoint paramPt(theLinPoint);
113  //tarjectory parameters
114  double field = theTrack.field()->inInverseGeV(thePredState.theState().position()).z();
115  double signTC = -theCharge;
116 
117  double thetaAtEP = thePredState.perigeeParameters().theta();
118  double phiAtEP = thePredState.perigeeParameters().phi();
119  double ptAtEP = thePredState.pt();
120  double transverseCurvatureAtEP = field / ptAtEP*signTC;
121 
122  double x_v = thePredState.theState().position().x();
123  double y_v = thePredState.theState().position().y();
124  double z_v = thePredState.theState().position().z();
125  double X = x_v - paramPt.x() - sin(phiAtEP) / transverseCurvatureAtEP;
126  double Y = y_v - paramPt.y() + cos(phiAtEP) / transverseCurvatureAtEP;
127  double SS = X*X + Y*Y;
128  double S = sqrt(SS);
129 
130  // The track parameters at the expansion point
131 
132  theExpandedParams[0] = transverseCurvatureAtEP;
133  theExpandedParams[1] = thetaAtEP;
134  theExpandedParams[3] = 1/transverseCurvatureAtEP - signTC * S;
135  double phiFEP;
136  if (std::abs(X)>std::abs(Y)) {
137  double signX = (X>0.0? +1.0:-1.0);
138  phiFEP = -signTC * signX*acos(signTC*Y/S);
139  } else {
140  phiFEP = asin(-signTC*X/S);
141  if ((signTC*Y)<0.0)
142  phiFEP = M_PI - phiFEP;
143  }
144  if (phiFEP>M_PI) phiFEP-= 2*M_PI;
145  theExpandedParams[2] = phiFEP;
146  theExpandedParams[4] = z_v - paramPt.z() -
147  (phiAtEP - theExpandedParams[2]) / tan(thetaAtEP)/transverseCurvatureAtEP;
148 
149  // The Jacobian: (all at the expansion point)
150  // [i,j]
151  // i = 0: rho , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
152  // j = 0: x_v, 1: y_v, 2: z_v
153 
154  thePositionJacobian(2,0) = - Y / (SS);
155  thePositionJacobian(2,1) = X / (SS);
156  thePositionJacobian(3,0) = - signTC * X / S;
157  thePositionJacobian(3,1) = - signTC * Y / S;
158  thePositionJacobian(4,0) = thePositionJacobian(2,0)/tan(thetaAtEP)/transverseCurvatureAtEP;
159  thePositionJacobian(4,1) = thePositionJacobian(2,1)/tan(thetaAtEP)/transverseCurvatureAtEP;
160  thePositionJacobian(4,2) = 1;
161 
162  // [i,j]
163  // i = 0: rho , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
164  // j = 0: rho, 1: theta, 2: phi_v
165  theMomentumJacobian(0,0) = 1;
166  theMomentumJacobian(1,1) = 1;
167 
168  theMomentumJacobian(2,0) = -
169  (X*cos(phiAtEP) + Y*sin(phiAtEP))/
170  (SS*transverseCurvatureAtEP*transverseCurvatureAtEP);
171 
172  theMomentumJacobian(2,2) = (Y*cos(phiAtEP) - X*sin(phiAtEP)) /
173  (SS*transverseCurvatureAtEP);
174 
175  theMomentumJacobian(3,0) =
176  (signTC * (Y*cos(phiAtEP) - X*sin(phiAtEP)) / S - 1)/
177  (transverseCurvatureAtEP*transverseCurvatureAtEP);
178 
179  theMomentumJacobian(3,2) = signTC *(X*cos(phiAtEP) + Y*sin(phiAtEP))/
180  (S*transverseCurvatureAtEP);
181 
182  theMomentumJacobian(4,0) = (phiAtEP - theExpandedParams[2]) /
183  tan(thetaAtEP)/(transverseCurvatureAtEP*transverseCurvatureAtEP)+
184  theMomentumJacobian(2,0) / tan(thetaAtEP)/transverseCurvatureAtEP;
185 
186  theMomentumJacobian(4,1) = (phiAtEP - theExpandedParams[2]) *
187  (1 + 1/(tan(thetaAtEP)*tan(thetaAtEP)))/transverseCurvatureAtEP;
188 
189  theMomentumJacobian(4,2) = (theMomentumJacobian(2,2) - 1) /
190  tan(thetaAtEP)/transverseCurvatureAtEP;
191 
192  // And finally the residuals:
193 
194  auto p = thePredState.theState().position();
195  AlgebraicVector3 expansionPoint(p.x(),p.y(),p.z());
196  AlgebraicVector3 momentumAtExpansionPoint( transverseCurvatureAtEP,thetaAtEP,phiAtEP);
197 
199  thePositionJacobian * expansionPoint -
200  theMomentumJacobian * momentumAtExpansionPoint );
201 
202 }
203 
204 
205 
206 
207 
209 {
210  GlobalPoint paramPt(theLinPoint);
211 
212  //tarjectory parameters
213  double thetaAtEP = thePredState.theState().momentum().theta();
214  double phiAtEP = thePredState.theState().momentum().phi();
215  double ptAtEP = thePredState.theState().momentum().perp();
216 
217  double x_v = thePredState.theState().position().x();
218  double y_v = thePredState.theState().position().y();
219  double z_v = thePredState.theState().position().z();
220  double X = x_v - paramPt.x();
221  double Y = y_v - paramPt.y();
222 
223  // The track parameters at the expansion point
224 
225  theExpandedParams(0) = 1 / ptAtEP;
226  theExpandedParams(1) = thetaAtEP;
227  theExpandedParams(2) = phiAtEP;
228  theExpandedParams(3) = X*sin(phiAtEP) - Y*cos(phiAtEP);
229  theExpandedParams(4) = z_v - paramPt.z() -
230  (X*cos(phiAtEP) + Y*sin(phiAtEP)) / tan(thetaAtEP);
231 
232  // The Jacobian: (all at the expansion point)
233  // [i,j]
234  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
235  // j = 0: x_v, 1: y_v, 2: z_v
236 
237  thePositionJacobian(3,0) = sin(phiAtEP);
238  thePositionJacobian(3,1) = - cos(phiAtEP);
239  thePositionJacobian(4,0) = - cos(phiAtEP)/tan(thetaAtEP);
240  thePositionJacobian(4,1) = - sin(phiAtEP)/tan(thetaAtEP);
241  thePositionJacobian(4,2) = 1;
242 
243  // [i,j]
244  // i = 0: rho = 1/pt , 1: theta, 2: phi_p, 3: epsilon, 4: z_p
245  // j = 0: rho = 1/pt , 1: theta, 2: phi_v
246 
247  theMomentumJacobian(0,0) = 1;
248  theMomentumJacobian(1,1) = 1;
249  theMomentumJacobian(2,2) = 1;
250 
251  theMomentumJacobian(3,2) = X*cos(phiAtEP) + Y*sin(phiAtEP);
252 
254  (1 + 1/(tan(thetaAtEP)*tan(thetaAtEP)));
255 
256  theMomentumJacobian(4,2) = (X*sin(phiAtEP) - Y*cos(phiAtEP))/tan(thetaAtEP);
257 
258  // And finally the residuals:
259 
260  auto p = thePredState.theState().position();
261  AlgebraicVector3 expansionPoint(p.x(),p.y(),p.z());
262  AlgebraicVector3 momentumAtExpansionPoint(1./ptAtEP,thetaAtEP,phiAtEP);
263 
265  thePositionJacobian * expansionPoint -
266  theMomentumJacobian * momentumAtExpansionPoint );
267 
268 }
269 
virtual AlgebraicVector3 predictedStateMomentumParameters() const
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const AlgebraicVector3 &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicSymMatrix66 &theCovarianceMatrix, const MagneticField *field)
T perp() const
Definition: PV3DBase.h:72
const AlgebraicMatrix53 & momentumJacobian() const
virtual RefCountedLinearizedTrackState stateWithNewLinearizationPoint(const GlobalPoint &newLP) const
ROOT::Math::SVector< double, N > AlgebraicVectorN
Common base class.
const FreeTrajectoryState & theState() const
virtual reco::TransientTrack track() const
virtual void checkParameters(AlgebraicVector5 &parameters) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
const MagneticField * field() const
bool operator==(LinearizedTrackState< 5 > &other) const
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
tuple result
Definition: mps_fire.py:84
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
#define unlikely(x)
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:39
TrajectoryStateClosestToPoint thePredState
const PerigeeTrajectoryParameters & perigeeParameters() const
T sqrt(T t)
Definition: SSEVec.h:18
const TrajectoryStateOnSurface theTSOS
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const AlgebraicMatrix53 & positionJacobian() const
ROOT::Math::SVector< double, 3 > AlgebraicVector3
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
GlobalVector momentum() const
#define M_PI
GlobalPoint position() const
ROOT::Math::SVector< double, 5 > AlgebraicVector5
virtual AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState &theRefittedState) const
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
T x() const
Definition: PV3DBase.h:62
virtual RefCountedRefittedTrackState createRefittedTrackState(const GlobalPoint &vertexPosition, const AlgebraicVector3 &vectorParameters, const AlgebraicSymMatrix66 &covarianceMatrix) const
virtual std::vector< ReferenceCountingPointer< LinearizedTrackState< 5 > > > components() const
const AlgebraicVector5 & constantTerm() const