CMS 3D CMS Logo

PerigeeRefittedTrackState.h
Go to the documentation of this file.
1 #ifndef PerigeeRefittedTrackState_H
2 #define PerigeeRefittedTrackState_H
3 
7 
15 class Surface;
16 class Propagator;
17 
19 
20 public:
21 
23 
25  const AlgebraicVector3 & aMomentumAtVertex,
26  const double aWeight = 1.) :
27  theState(tscp), momentumAtVertex(aMomentumAtVertex), theWeight(aWeight) {}
28 
30 
36  {return theState.theState();}
37 
42  const Surface & surface) const;
43 
49  const Surface & surface, const Propagator & propagator) const;
50 
57  virtual AlgebraicVector5 parameters() const
58  {return theState.perigeeParameters().vector();}
59 
66 
71  virtual GlobalPoint position() const
72  {return theState.referencePoint();}
73 
79  virtual AlgebraicVector3 momentumVector() const;
80 
84  virtual double weight() const {return theWeight;}
85 
91  (const double newWeight) const;
92 
93  virtual std::vector<ReferenceCountingPointer<RefittedTrackState<5> > > components() const;
94 
95  virtual reco::TransientTrack transientTrack() const;
96 
97 private:
98 
101  double theWeight;
102 };
103 #endif
const AlgebraicVector5 & vector() const
virtual FreeTrajectoryState freeTrajectoryState() const
const PerigeeTrajectoryError & perigeeError() const
const FreeTrajectoryState & theState() const
virtual double weight() const
virtual AlgebraicSymMatrix55 covariance() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
ReferenceCountingPointer< RefittedTrackState< 5 > > RefCountedRefittedTrackState
virtual AlgebraicVector5 parameters() const
virtual ReferenceCountingPointer< RefittedTrackState< 5 > > stateWithNewWeight(const double newWeight) const
virtual AlgebraicVector3 momentumVector() const
const PerigeeTrajectoryParameters & perigeeParameters() const
ROOT::Math::SVector< double, 3 > AlgebraicVector3
virtual std::vector< ReferenceCountingPointer< RefittedTrackState< 5 > > > components() const
PerigeeRefittedTrackState(const TrajectoryStateClosestToPoint &tscp, const AlgebraicVector3 &aMomentumAtVertex, const double aWeight=1.)
const AlgebraicSymMatrix55 & covarianceMatrix() const
TrajectoryStateClosestToPoint theState
ROOT::Math::SVector< double, 5 > AlgebraicVector5
virtual TrajectoryStateOnSurface trajectoryStateOnSurface(const Surface &surface) const
const GlobalPoint & referencePoint() const
virtual reco::TransientTrack transientTrack() const
virtual GlobalPoint position() const