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 const GlobalTrajectoryParameters& parameters() const {
00100 return theGlobalParameters;
00101 }
00102 const CartesianTrajectoryError& cartesianError() const {
00103 if (!hasError()) throw TrajectoryStateException(
00104 "FreeTrajectoryState: attempt to access errors when none available");
00105 if (!theCartesianErrorValid)
00106 createCartesianError();
00107 return theCartesianError;
00108 }
00109 const CurvilinearTrajectoryError& curvilinearError() const {
00110 if (!hasError()) throw TrajectoryStateException(
00111 "FreeTrajectoryState: attempt to access errors when none available");
00112 if (!theCurvilinearErrorValid)
00113 createCurvilinearError();
00114 return theCurvilinearError;
00115 }
00116 void rescaleError(double factor) {
00117 bool zeroField = parameters().magneticFieldInInverseGeV(position()).mag()==0;
00118 if (zeroField) {
00119 if (theCartesianErrorValid){
00120 if (!theCurvilinearErrorValid) createCurvilinearError();
00121 theCurvilinearError.zeroFieldScaling(factor*factor);
00122 createCartesianError();
00123 }else
00124 if (theCurvilinearErrorValid) theCurvilinearError.zeroFieldScaling(factor*factor);
00125 } else{
00126 if (theCartesianErrorValid){
00127 theCartesianError *= (factor*factor);
00128 }
00129 if (theCurvilinearErrorValid){
00130 theCurvilinearError *= (factor*factor);
00131 }
00132 }
00133 }
00134
00135 void setCartesianError(const CartesianTrajectoryError &err) {
00136 theCartesianError = err; theCartesianErrorValid = true;
00137 }
00138 void setCartesianError(const AlgebraicSymMatrix66 &err) {
00139 theCartesianError = CartesianTrajectoryError(err); theCartesianErrorValid = true;
00140 }
00141 void setCurvilinearError(const CurvilinearTrajectoryError &err) {
00142 theCurvilinearError = err; theCurvilinearErrorValid = true;
00143 }
00144 void setCurvilinearError(const AlgebraicSymMatrix55 &err) {
00145 theCurvilinearError = CurvilinearTrajectoryError(err); theCurvilinearErrorValid = true;
00146 }
00147
00148 bool canReach(double radius) const;
00149 private:
00150
00151 void createCartesianError() const;
00152
00153 void createCurvilinearError() const;
00154
00155 private:
00156
00157 GlobalTrajectoryParameters theGlobalParameters;
00158 CartesianTrajectoryError theCartesianError;
00159 CurvilinearTrajectoryError theCurvilinearError;
00160 bool theCartesianErrorValid;
00161 bool theCurvilinearErrorValid;
00162
00163 };
00164
00165 std::ostream& operator<<(std::ostream& os, const FreeTrajectoryState& fts);
00166
00167 #endif