CMS 3D CMS Logo

PerigeeLinearizedTrackState.h
Go to the documentation of this file.
1 #ifndef PerigeeLinearizedTrackState_H
2 #define PerigeeLinearizedTrackState_H
3 
8 #include "Math/SMatrix.h"
10 
12 
35 public:
41 
46  RefCountedLinearizedTrackState stateWithNewLinearizationPoint(const GlobalPoint& newLP) const override;
47 
51  const GlobalPoint& linearizationPoint() const override { return theLinPoint; }
52 
53  reco::TransientTrack track() const override { return theTrack; }
54 
55  const TrajectoryStateOnSurface state() const { return theTSOS; }
56 
60  const AlgebraicVector5& constantTerm() const override;
61 
65  const AlgebraicMatrix53& positionJacobian() const override;
66 
70  const AlgebraicMatrix53& momentumJacobian() const override;
71 
74  const AlgebraicVector5& parametersFromExpansion() const override;
75 
81 
86 
91 
96  AlgebraicSymMatrix55 predictedStateWeight(int& error) const override;
97 
101  AlgebraicSymMatrix55 predictedStateError() const override;
102 
107 
108  // /** Method returning the impact point measurement
109  // */
110  // ImpactPointMeasurement impactPointMeasurement() const;
111 
112  TrackCharge charge() const override { return theCharge; }
113 
114  bool hasError() const override;
115 
116  bool operator==(LinearizedTrackState<5>& other) const override;
117 
119 
124  const AlgebraicVector3& vectorParameters,
125  const AlgebraicSymMatrix66& covarianceMatrix) const override;
126 
127  AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState& theRefittedState) const override;
128 
129  double weightInMixture() const override { return theTSOS.weight(); }
130 
131  void inline checkParameters(AlgebraicVector5& parameters) const override;
132 
133  std::vector<ReferenceCountingPointer<LinearizedTrackState<5> > > components() const override;
134 
135  bool isValid() const override;
136 
137 private:
143  const TrajectoryStateOnSurface& tsos)
145 
148  void computeJacobians() const;
151  void computeChargedJacobians() const;
154  void computeNeutralJacobians() const;
157  void compute3DImpactPoint() const;
158 
162 
167 
169 
171  mutable bool jacobiansAvailable;
172 };
173 
178  if
180  return theConstantTerm;
181 }
182 
187  if
189  return thePositionJacobian;
190 }
191 
196  if
198  return theMomentumJacobian;
199 }
200 
204  if
206  return theExpandedParams;
207 }
208 
214  if
216  return thePredState;
217 }
218 
220  if
222  return thePredState.hasError();
223 }
224 
226  if (!jacobiansAvailable)
229 }
230 
232  if
235  return AlgebraicVector3(v[0], v[1], v[2]);
236 }
237 
239  if
241  if (!thePredState.isValid()) {
242  error = 1;
243  return AlgebraicSymMatrix55();
244  }
245  return thePredState.perigeeError().weightMatrix(error);
246 }
247 
249  if
252 }
253 
255  if
258 }
259 
261  if
262  UNLIKELY(!theTSOS.isValid()) return false;
263  if
265  return jacobiansAvailable;
266 }
267 
268 #endif
RefCountedRefittedTrackState createRefittedTrackState(const GlobalPoint &vertexPosition, const AlgebraicVector3 &vectorParameters, const AlgebraicSymMatrix66 &covarianceMatrix) const override
const AlgebraicVector5 & vector() const
AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState &theRefittedState) const override
const PerigeeTrajectoryError & perigeeError() const
void compute3DImpactPoint() const
const GlobalPoint & linearizationPoint() const override
RefCountedLinearizedTrackState stateWithNewLinearizationPoint(const GlobalPoint &newLP) const override
const TrajectoryStateOnSurface state() const
AlgebraicVector3 predictedStateMomentumParameters() const override
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
int TrackCharge
Definition: TrackCharge.h:4
TrajectoryStateClosestToPoint thePredState
bool operator==(LinearizedTrackState< 5 > &other) const override
const PerigeeTrajectoryParameters & perigeeParameters() const
ROOT::Math::SMatrix< double, 5, 3, ROOT::Math::MatRepStd< double, 5, 3 > > AlgebraicMatrix53
const TrajectoryStateClosestToPoint & predictedState() const
const TrajectoryStateOnSurface theTSOS
const AlgebraicVector5 & parametersFromExpansion() const override
AlgebraicSymMatrix33 predictedStateMomentumError() const override
ROOT::Math::SVector< double, 5 > AlgebraicVector5
std::vector< ReferenceCountingPointer< LinearizedTrackState< 5 > > > components() const override
const AlgebraicSymMatrix55 & covarianceMatrix() const
PerigeeLinearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track, const TrajectoryStateOnSurface &tsos)
const AlgebraicSymMatrix55 & weightMatrix(int &error) const
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66
reco::TransientTrack track() const override
AlgebraicSymMatrix55 predictedStateError() const override
const AlgebraicMatrix53 & momentumJacobian() const override
const AlgebraicVector5 & constantTerm() const override
const AlgebraicMatrix53 & positionJacobian() const override
TrackCharge charge() const override
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
#define UNLIKELY(x)
Definition: Likely.h:21
AlgebraicSymMatrix55 predictedStateWeight(int &error) const override
void checkParameters(AlgebraicVector5 &parameters) const override
AlgebraicVector5 predictedStateParameters() const override