CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FreeTrajectoryState.h
Go to the documentation of this file.
1 #ifndef _TRACKER_FREETRAJECTORYSTATE_H_
2 #define _TRACKER_FREETRAJECTORYSTATE_H_
3 
4 // base trajectory state class
5 // track parameters and error covariance matrix
6 
12 
13 #include <iosfwd>
14 
16 
17 
20 
21 
31 class FreeTrajectoryState : public BlockWipedAllocated<FreeTrajectoryState> {
32 public:
33 // construst
34 //default constructor - needed for Persistency
35 
38 
39  FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters) :
40  theGlobalParameters(aGlobalParameters),
42  {}
43 
45  const GlobalVector& aP,
46  TrackCharge aCharge,
47  const MagneticField* fieldProvider) :
48  theGlobalParameters(aX, aP, aCharge, fieldProvider),
50  {}
51 
52 
54  const CurvilinearTrajectoryError& aCurvilinearError) :
55  theGlobalParameters(aGlobalParameters),
56  theCurvilinearError(aCurvilinearError)
57  {}
58 
59 
60 
62  const CartesianTrajectoryError& aCartesianError) :
63  theGlobalParameters(aGlobalParameters)
64  { createCurvilinearError(aCartesianError); }
65 
68  const CurvilinearTrajectoryError& aCurvilinearError) :
69  theGlobalParameters(aGlobalParameters),
70  theCurvilinearError(aCurvilinearError)
71  {}
72 
73 // access
74 // propagate access to parameters
77  }
80  }
81  TrackCharge charge() const {
82  return theGlobalParameters.charge();
83  }
84  double signedInverseMomentum() const {
86  }
87  double transverseCurvature() const {
89  }
90 
91 // direct access
92 
94 
95  bool hasError() const {
96  return hasCurvilinearError();
97  }
98 
99 
101  return theGlobalParameters;
102  }
103 
104 
106  if unlikely(!hasError()) missingError();
107  CartesianTrajectoryError aCartesianError;
108  createCartesianError(aCartesianError);
109  return aCartesianError;
110  }
111 
113  if unlikely(!hasError()) missingError();
114  return theCurvilinearError;
115  }
116 
117  void rescaleError(double factor);
118 
119 
122  }
125  }
126 
128  theCurvilinearError = err;
129  }
132  }
133 
134 
135 // properties
136  bool canReach(double radius) const;
137 private:
138 
139 
140  void missingError() const; // dso_internal;
141 
142 // convert curvilinear errors to cartesian
143  void createCartesianError(CartesianTrajectoryError & aCartesianError) const; // dso_internal;
144 
145 
146 // convert cartesian errors to curvilinear
147  void createCurvilinearError(CartesianTrajectoryError const & aCartesianError) const; // dso_internal;
148 
149 private:
150 
153 
154 };
155 
156 std::ostream& operator<<(std::ostream& os, const FreeTrajectoryState& fts);
157 
158 #endif
CartesianTrajectoryError cartesianError() const
const GlobalTrajectoryParameters & parameters() const
bool hasCurvilinearError() const
void setCartesianError(const AlgebraicSymMatrix66 &err)
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66
CurvilinearTrajectoryError theCurvilinearError
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TrackCharge charge() const
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
const CurvilinearTrajectoryError & curvilinearError() const
void setCurvilinearError(const AlgebraicSymMatrix55 &err)
FreeTrajectoryState(const GlobalTrajectoryParameters &aGlobalParameters, const CartesianTrajectoryError &, const CurvilinearTrajectoryError &aCurvilinearError)
FreeTrajectoryState(const GlobalTrajectoryParameters &aGlobalParameters, const CartesianTrajectoryError &aCartesianError)
#define unlikely(x)
Definition: Likely.h:21
FreeTrajectoryState(const GlobalPoint &aX, const GlobalVector &aP, TrackCharge aCharge, const MagneticField *fieldProvider)
GlobalTrajectoryParameters theGlobalParameters
int TrackCharge
Definition: TrackCharge.h:4
void createCurvilinearError(CartesianTrajectoryError const &aCartesianError) const
FreeTrajectoryState(const GlobalTrajectoryParameters &aGlobalParameters, const CurvilinearTrajectoryError &aCurvilinearError)
GlobalVector momentum() const
void setCurvilinearError(const CurvilinearTrajectoryError &err)
void setCartesianError(const CartesianTrajectoryError &err)
GlobalPoint position() const
void rescaleError(double factor)
void createCartesianError(CartesianTrajectoryError &aCartesianError) const
bool canReach(double radius) const
FreeTrajectoryState(const GlobalTrajectoryParameters &aGlobalParameters)
double transverseCurvature() const
double signedInverseMomentum() const