CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
13 
35 class PerigeeLinearizedTrackState GCC11_FINAL : public LinearizedTrackState<5> {
36 
37 
38 public:
39 
45 
51  (const GlobalPoint & newLP) const;
52 
53 
57  const GlobalPoint & linearizationPoint() const { return theLinPoint; }
58 
59  virtual reco::TransientTrack track() const { return theTrack; }
60 
61  const TrajectoryStateOnSurface state() const { return theTSOS; }
62 
66  const AlgebraicVector5 & constantTerm() const;
67 
71  const AlgebraicMatrix53 & positionJacobian() const;
72 
76  const AlgebraicMatrix53 & momentumJacobian() const;
77 
81 
86  const TrajectoryStateClosestToPoint & predictedState() const;
87 
92 
97 
103 
108 
113 
114 // /** Method returning the impact point measurement
115 // */
116 // ImpactPointMeasurement impactPointMeasurement() const;
117 
118  TrackCharge charge() const {return theCharge;}
119 
120  bool hasError() const;
121 
122  bool operator ==(LinearizedTrackState<5> & other)const;
123 
125 
129  virtual RefCountedRefittedTrackState createRefittedTrackState(
130  const GlobalPoint & vertexPosition,
131  const AlgebraicVector3 & vectorParameters,
132  const AlgebraicSymMatrix66 & covarianceMatrix) const;
133 
134 
136  const RefCountedRefittedTrackState & theRefittedState) const;
137 
138  virtual double weightInMixture() const {return theTSOS.weight();}
139 
140  virtual void inline checkParameters(AlgebraicVector5 & parameters) const;
141 
142  virtual std::vector<ReferenceCountingPointer<LinearizedTrackState<5> > > components() const;
143 
144  virtual bool isValid() const;
145 
146 private:
147 
152  const TrajectoryStateOnSurface& tsos)
153  : theLinPoint(linP), theTrack(track),
154  theTSOS(tsos), theCharge(theTrack.charge()), jacobiansAvailable(false) {}
155 
158  void computeJacobians() const;
161  void computeChargedJacobians() const;
164  void computeNeutralJacobians() const;
167  void compute3DImpactPoint() const;
168 
172 
173  mutable AlgebraicMatrix53 thePositionJacobian, theMomentumJacobian;
177 
178  TSCPBuilderNoMaterial builder;
179 
181  mutable bool jacobiansAvailable;
182 
183 };
184 
185 
186 
187 
191 inline
192 const AlgebraicVector5 & PerigeeLinearizedTrackState::constantTerm() const
193 {
194  if unlikely(!jacobiansAvailable) computeJacobians();
195  return theConstantTerm;
196 }
197 
201 inline
202 const AlgebraicMatrix53 & PerigeeLinearizedTrackState::positionJacobian() const
203 {
204  if unlikely(!jacobiansAvailable) computeJacobians();
205  return thePositionJacobian;
206 }
207 
211 inline
212 const AlgebraicMatrix53 & PerigeeLinearizedTrackState::momentumJacobian() const
213 {
214  if unlikely(!jacobiansAvailable) computeJacobians();
215  return theMomentumJacobian;
216 }
217 
220 inline
221 const AlgebraicVector5 & PerigeeLinearizedTrackState::parametersFromExpansion() const
222 {
223  if unlikely(!jacobiansAvailable) computeJacobians();
224  return theExpandedParams;
225 }
226 
231 inline
232 const TrajectoryStateClosestToPoint & PerigeeLinearizedTrackState::predictedState() const
233 {
234  if unlikely(!jacobiansAvailable) computeJacobians();
235  return thePredState;
236 }
237 
238 
239 inline
240 bool PerigeeLinearizedTrackState::hasError() const
241 {
242  if unlikely(!jacobiansAvailable) computeJacobians();
243  return thePredState.hasError();
244 }
245 
246 inline
247 AlgebraicVector5 PerigeeLinearizedTrackState::predictedStateParameters() const
248 {
249  if (!jacobiansAvailable) computeJacobians();
250  return thePredState.perigeeParameters().vector();
251 }
252 
253 inline
254 AlgebraicVector3 PerigeeLinearizedTrackState::predictedStateMomentumParameters() const
255 {
256  if unlikely(!jacobiansAvailable) computeJacobians();
257  auto v= thePredState.perigeeParameters().vector();
258  return AlgebraicVector3(v[0],v[1],v[2]);
259 
260 }
261 
262 inline
263 AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateWeight(int & error) const
264 {
265  if unlikely(!jacobiansAvailable) computeJacobians();
266  if (!thePredState.isValid()) {
267  error = 1;
268  return AlgebraicSymMatrix55();
269  }
270  return thePredState.perigeeError().weightMatrix(error);
271 }
272 
273 inline
274 AlgebraicSymMatrix55 PerigeeLinearizedTrackState::predictedStateError() const
275 {
276  if unlikely(!jacobiansAvailable) computeJacobians();
277  return thePredState.perigeeError().covarianceMatrix();
278 }
279 
280 inline
281 AlgebraicSymMatrix33 PerigeeLinearizedTrackState::predictedStateMomentumError() const
282 {
283  if unlikely(!jacobiansAvailable) computeJacobians();
284  return thePredState.perigeeError().covarianceMatrix().Sub<AlgebraicSymMatrix33>(0,2);
285 }
286 
287 inline
288 bool PerigeeLinearizedTrackState::isValid() const
289 {
290  if unlikely(!theTSOS.isValid()) return false;
291  if unlikely(!jacobiansAvailable) computeJacobians();
292  return jacobiansAvailable;
293 }
294 
295 
296 
297 #endif
virtual bool operator==(LinearizedTrackState< N > &other) const =0
dictionary parameters
Definition: Parameters.py:2
AlgebraicMatrix53 thePositionJacobian
virtual const AlgebraicVectorN & constantTerm() const =0
AlgebraicVector5 theExpandedParams
virtual const AlgebraicMatrixNM & momentumJacobian() const =0
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
virtual ReferenceCountingPointer< LinearizedTrackState< N > > stateWithNewLinearizationPoint(const GlobalPoint &newLP) const =0
GloballyPositioned< float >::GlobalPoint GlobalPoint
Definition: MagVolume.h:20
const TrajectoryStateOnSurface theTSOS
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66
PerigeeLinearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track, const TrajectoryStateOnSurface &tsos)
const GlobalPoint & linearizationPoint() const
TrajectoryStateClosestToPoint thePredState
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
double charge(const std::vector< uint8_t > &Ampls)
reco::TransientTrack theTrack
const TrajectoryStateOnSurface state() const
virtual AlgebraicVectorM predictedStateMomentumParameters() const =0
#define unlikely(x)
Definition: Likely.h:21
virtual const AlgebraicVectorN & parametersFromExpansion() const =0
int TrackCharge
Definition: TrackCharge.h:4
virtual RefCountedRefittedTrackState createRefittedTrackState(const GlobalPoint &vertexPosition, const AlgebraicVectorM &vectorParameters, const AlgebraicSymMatrixOO &covarianceMatrix) const =0
ROOT::Math::SMatrix< double, 5, 3, ROOT::Math::MatRepStd< double, 5, 3 > > AlgebraicMatrix53
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
virtual AlgebraicSymMatrixNN predictedStateError() const =0
virtual AlgebraicVectorN predictedStateParameters() const =0
virtual AlgebraicVectorN refittedParamFromEquation(const RefCountedRefittedTrackState &theRefittedState) const =0
ROOT::Math::SVector< double, 3 > AlgebraicVector3
virtual reco::TransientTrack track() const
virtual void checkParameters(AlgebraicVectorN &parameters) const
virtual bool isValid() const
virtual double weightInMixture() const
ROOT::Math::SVector< double, 5 > AlgebraicVector5
string const
Definition: compareJSON.py:14
virtual bool hasError() const =0
TrackCharge charge() const
return(e1-e2)*(e1-e2)+dp *dp
virtual std::vector< ReferenceCountingPointer< LinearizedTrackState< N > > > components() const =0
TSCPBuilderNoMaterial builder
if(dp >Float(M_PI)) dp-
virtual AlgebraicSymMatrixNN predictedStateWeight(int &error) const =0
virtual const AlgebraicMatrixN3 & positionJacobian() const =0
virtual AlgebraicSymMatrixMM predictedStateMomentumError() const =0
volatile std::atomic< bool > shutdown_flag false
AlgebraicVector5 theConstantTerm
Unlimited (trivial) bounds.