Go to the documentation of this file.00001 #ifndef _TRACKER_FREETRAJECTORYSTATE_H_
00002 #define _TRACKER_FREETRAJECTORYSTATE_H_
00003
00004
00005
00006
00007 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00008 #include "TrackingTools/TrajectoryParametrization/interface/CartesianTrajectoryError.h"
00009 #include "TrackingTools/TrajectoryParametrization/interface/CurvilinearTrajectoryError.h"
00010 #include "DataFormats/TrajectoryState/interface/TrackCharge.h"
00011 #include "TrackingTools/TrajectoryParametrization/interface/TrajectoryStateExceptions.h"
00012
00013 #include <iosfwd>
00014
00024 class FreeTrajectoryState {
00025 public:
00026
00027
00028
00029 FreeTrajectoryState():
00030 theCartesianErrorValid(false),
00031 theCurvilinearErrorValid(false) {};
00032
00033 FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters) :
00034 theGlobalParameters(aGlobalParameters),
00035 theCartesianErrorValid(false),
00036 theCurvilinearErrorValid(false)
00037 {
00038 }
00039
00040 FreeTrajectoryState(const GlobalPoint& aX,
00041 const GlobalVector& aP,
00042 TrackCharge aCharge,
00043 const MagneticField* fieldProvider) :
00044 theGlobalParameters(aX, aP, aCharge, fieldProvider),
00045 theCartesianErrorValid(false),
00046 theCurvilinearErrorValid(false)
00047 {
00048 }
00049
00050 FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00051 const CartesianTrajectoryError& aCartesianError) :
00052 theGlobalParameters(aGlobalParameters),
00053 theCartesianError(aCartesianError),
00054 theCartesianErrorValid(true),
00055 theCurvilinearErrorValid(false)
00056 {
00057 }
00058 FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00059 const CurvilinearTrajectoryError& aCurvilinearError) :
00060 theGlobalParameters(aGlobalParameters),
00061 theCurvilinearError(aCurvilinearError),
00062 theCartesianErrorValid(false),
00063 theCurvilinearErrorValid(true)
00064 {
00065 }
00066 FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00067 const CartesianTrajectoryError& aCartesianError,
00068 const CurvilinearTrajectoryError& aCurvilinearError) :
00069 theGlobalParameters(aGlobalParameters),
00070 theCartesianError(aCartesianError),
00071 theCurvilinearError(aCurvilinearError),
00072 theCartesianErrorValid(true),
00073 theCurvilinearErrorValid(true)
00074 {
00075 }
00076
00077
00078 GlobalPoint position() const {
00079 return theGlobalParameters.position();
00080 }
00081 GlobalVector momentum() const {
00082 return theGlobalParameters.momentum();
00083 }
00084 TrackCharge charge() const {
00085 return theGlobalParameters.charge();
00086 }
00087 double signedInverseMomentum() const {
00088 return theGlobalParameters.signedInverseMomentum();
00089 }
00090 double transverseCurvature() const {
00091 return theGlobalParameters.transverseCurvature();
00092 }
00093
00094 bool hasCartesianError() const {return theCartesianErrorValid;}
00095 bool hasCurvilinearError() const {return theCurvilinearErrorValid;}
00096 bool hasError() const {
00097 return theCurvilinearErrorValid || theCartesianErrorValid;
00098 }
00099 static void missingError();
00100
00101 const GlobalTrajectoryParameters& parameters() const {
00102 return theGlobalParameters;
00103 }
00104 const CartesianTrajectoryError& cartesianError() const {
00105 if (!hasError()) missingError();
00106 if (!theCartesianErrorValid)
00107 createCartesianError();
00108 return theCartesianError;
00109 }
00110 const CurvilinearTrajectoryError& curvilinearError() const {
00111 if (!hasError()) missingError();
00112 if (!theCurvilinearErrorValid)
00113 createCurvilinearError();
00114 return theCurvilinearError;
00115 }
00116
00117 void rescaleError(double factor);
00118
00119 void setCartesianError(const CartesianTrajectoryError &err) {
00120 theCartesianError = err; theCartesianErrorValid = true;
00121 }
00122 void setCartesianError(const AlgebraicSymMatrix66 &err) {
00123 theCartesianError = CartesianTrajectoryError(err); theCartesianErrorValid = true;
00124 }
00125 void setCurvilinearError(const CurvilinearTrajectoryError &err) {
00126 theCurvilinearError = err; theCurvilinearErrorValid = true;
00127 }
00128 void setCurvilinearError(const AlgebraicSymMatrix55 &err) {
00129 theCurvilinearError = CurvilinearTrajectoryError(err); theCurvilinearErrorValid = true;
00130 }
00131
00132 bool canReach(double radius) const;
00133 private:
00134
00135 void createCartesianError() const;
00136
00137 void createCurvilinearError() const;
00138
00139 private:
00140
00141 GlobalTrajectoryParameters theGlobalParameters;
00142 CartesianTrajectoryError theCartesianError;
00143 CurvilinearTrajectoryError theCurvilinearError;
00144 bool theCartesianErrorValid;
00145 bool theCurvilinearErrorValid;
00146
00147 };
00148
00149 std::ostream& operator<<(std::ostream& os, const FreeTrajectoryState& fts);
00150
00151 #endif