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"
11 
34 public:
40 
46 
50  const GlobalPoint& linearizationPoint() const override { return theLinPoint; }
51 
52  reco::TransientTrack track() const override { return theTrack; }
53 
54  const TrajectoryStateOnSurface state() const { return theTSOS; }
55 
59  const AlgebraicVector5& constantTerm() const override;
60 
64  const AlgebraicMatrix53& positionJacobian() const override;
65 
69  const AlgebraicMatrix53& momentumJacobian() const override;
70 
73  const AlgebraicVector5& parametersFromExpansion() const override;
74 
80 
85 
90 
95  AlgebraicSymMatrix55 predictedStateWeight(int& error) const override;
96 
100  AlgebraicSymMatrix55 predictedStateError() const override;
101 
106 
107  // /** Method returning the impact point measurement
108  // */
109  // ImpactPointMeasurement impactPointMeasurement() const;
110 
111  TrackCharge charge() const override { return theCharge; }
112 
113  bool hasError() const override;
114 
115  bool operator==(LinearizedTrackState<5>& other) const override;
116 
118 
123  const AlgebraicVector3& vectorParameters,
124  const AlgebraicSymMatrix66& covarianceMatrix) const override;
125 
126  AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState& theRefittedState) const override;
127 
128  double weightInMixture() const override { return theTSOS.weight(); }
129 
130  void inline checkParameters(AlgebraicVector5& parameters) const override;
131 
132  std::vector<ReferenceCountingPointer<LinearizedTrackState<5> > > components() const override;
133 
134  bool isValid() const override;
135 
136 private:
142  const TrajectoryStateOnSurface& tsos)
144 
147  void computeJacobians() const;
150  void computeChargedJacobians() const;
153  void computeNeutralJacobians() const;
156  void compute3DImpactPoint() const;
157 
161 
166 
168 
170  mutable bool jacobiansAvailable;
171 };
172 
179  return theConstantTerm;
180 }
181 
188  return thePositionJacobian;
189 }
190 
197  return theMomentumJacobian;
198 }
199 
205  return theExpandedParams;
206 }
207 
215  return thePredState;
216 }
217 
221  return thePredState.hasError();
222 }
223 
225  if (!jacobiansAvailable)
228 }
229 
234  return AlgebraicVector3(v[0], v[1], v[2]);
235 }
236 
240  if (!thePredState.isValid()) {
241  error = 1;
242  return AlgebraicSymMatrix55();
243  }
245 }
246 
251 }
252 
257 }
258 
260  if UNLIKELY (!theTSOS.isValid())
261  return false;
264  return jacobiansAvailable;
265 }
266 
267 #endif
std::vector< ReferenceCountingPointer< LinearizedTrackState< 5 > > > components() const override
const AlgebraicVector5 & vector() const
TrackCharge charge() const override
const AlgebraicSymMatrix55 & covarianceMatrix() const
void checkParameters(AlgebraicVector5 &parameters) const override
const AlgebraicVector5 & parametersFromExpansion() const override
void compute3DImpactPoint() const
AlgebraicSymMatrix55 predictedStateWeight(int &error) const override
const TrajectoryStateClosestToPoint & predictedState() const
AlgebraicVector3 predictedStateMomentumParameters() const override
const TrajectoryStateOnSurface state() const
RefCountedLinearizedTrackState stateWithNewLinearizationPoint(const GlobalPoint &newLP) const override
const PerigeeTrajectoryError & perigeeError() const
const AlgebraicVector5 & constantTerm() const override
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
int TrackCharge
Definition: TrackCharge.h:4
TrajectoryStateClosestToPoint thePredState
ROOT::Math::SMatrix< double, 5, 3, ROOT::Math::MatRepStd< double, 5, 3 > > AlgebraicMatrix53
bool operator==(LinearizedTrackState< 5 > &other) const override
const TrajectoryStateOnSurface theTSOS
ROOT::Math::SVector< double, 5 > AlgebraicVector5
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
const PerigeeTrajectoryParameters & perigeeParameters() const
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
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
reco::TransientTrack track() const override
AlgebraicSymMatrix55 predictedStateError() const override
AlgebraicVector5 refittedParamFromEquation(const RefCountedRefittedTrackState &theRefittedState) const override
ROOT::Math::SVector< double, 3 > AlgebraicVector3
AlgebraicVector5 predictedStateParameters() const override
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
AlgebraicSymMatrix33 predictedStateMomentumError() const override
#define UNLIKELY(x)
Definition: Likely.h:21
const AlgebraicSymMatrix55 & weightMatrix(int &error) const
RefCountedRefittedTrackState createRefittedTrackState(const GlobalPoint &vertexPosition, const AlgebraicVector3 &vectorParameters, const AlgebraicSymMatrix66 &covarianceMatrix) const override
const GlobalPoint & linearizationPoint() const override