CMS 3D CMS Logo

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