CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h

Go to the documentation of this file.
00001 #ifndef _TRACKER_FREETRAJECTORYSTATE_H_
00002 #define _TRACKER_FREETRAJECTORYSTATE_H_
00003 
00004 // base trajectory state class
00005 // track parameters and error covariance matrix
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 
00015 #include "DataFormats/GeometrySurface/interface/BlockWipedAllocator.h"
00016 
00017 
00018 #include "FWCore/Utilities/interface/Visibility.h"
00019 #include "FWCore/Utilities/interface/Likely.h"
00020 
00021 
00031 class FreeTrajectoryState : public BlockWipedAllocated<FreeTrajectoryState> {
00032 public:
00033 // construst
00034 //default constructor - needed for Persistency
00035 
00036   FreeTrajectoryState():  
00037     theCurvilinearError(InvalidError()) {}
00038   
00039   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters) :
00040   theGlobalParameters(aGlobalParameters),  
00041   theCurvilinearError(InvalidError())
00042   {}
00043   
00044   FreeTrajectoryState(const GlobalPoint& aX,
00045                       const GlobalVector& aP,
00046                       TrackCharge aCharge, 
00047                       const MagneticField* fieldProvider) :
00048     theGlobalParameters(aX, aP, aCharge, fieldProvider),  
00049     theCurvilinearError(InvalidError())
00050   {}
00051   
00052   
00053   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00054                       const CurvilinearTrajectoryError& aCurvilinearError) :
00055     theGlobalParameters(aGlobalParameters),
00056     theCurvilinearError(aCurvilinearError)
00057   {}
00058   
00059   
00060   
00061   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00062                       const CartesianTrajectoryError& aCartesianError) :
00063     theGlobalParameters(aGlobalParameters)
00064   { createCurvilinearError(aCartesianError);  }
00065 
00066   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00067                       const CartesianTrajectoryError&,
00068                       const CurvilinearTrajectoryError& aCurvilinearError) :
00069     theGlobalParameters(aGlobalParameters),
00070     theCurvilinearError(aCurvilinearError)  
00071   {}
00072   
00073 // access
00074 // propagate access to parameters
00075   GlobalPoint position() const {
00076     return theGlobalParameters.position();
00077   }
00078   GlobalVector momentum() const {
00079     return theGlobalParameters.momentum();
00080   }
00081   TrackCharge charge() const {
00082     return theGlobalParameters.charge();
00083   }
00084   double signedInverseMomentum() const {
00085     return theGlobalParameters.signedInverseMomentum();
00086   }
00087   double transverseCurvature() const {
00088     return theGlobalParameters.transverseCurvature();
00089   }
00090 
00091 // direct access
00092 
00093   bool hasCurvilinearError() const {return theCurvilinearError.valid();}
00094 
00095   bool hasError() const {
00096     return hasCurvilinearError();
00097   }
00098 
00099 
00100   const GlobalTrajectoryParameters& parameters() const {
00101     return theGlobalParameters;
00102   }
00103 
00104 
00105   CartesianTrajectoryError cartesianError() const {
00106     if unlikely(!hasError()) missingError();
00107     CartesianTrajectoryError aCartesianError;
00108     createCartesianError(aCartesianError);
00109     return aCartesianError;
00110   }
00111 
00112   const CurvilinearTrajectoryError& curvilinearError() const {
00113     if  unlikely(!hasError()) missingError();
00114     return theCurvilinearError;
00115   }
00116 
00117   void rescaleError(double factor);
00118 
00119 
00120   void setCartesianError(const CartesianTrajectoryError &err) {
00121     createCurvilinearError(err);
00122   }
00123   void setCartesianError(const AlgebraicSymMatrix66 &err) {
00124     createCurvilinearError(CartesianTrajectoryError(err)); 
00125   }
00126 
00127   void setCurvilinearError(const CurvilinearTrajectoryError &err) {
00128         theCurvilinearError = err;
00129   }
00130   void setCurvilinearError(const AlgebraicSymMatrix55 &err) {
00131         theCurvilinearError = CurvilinearTrajectoryError(err); 
00132   }
00133 
00134 
00135 // properties
00136   bool canReach(double radius) const;
00137 private:
00138 
00139 
00140   void missingError() const; // dso_internal;
00141 
00142 // convert curvilinear errors to cartesian
00143   void createCartesianError(CartesianTrajectoryError & aCartesianError) const; // dso_internal;
00144 
00145 
00146 // convert cartesian errors to curvilinear
00147   void createCurvilinearError(CartesianTrajectoryError const & aCartesianError) const; // dso_internal;
00148 
00149 private:
00150 
00151   GlobalTrajectoryParameters  theGlobalParameters;
00152   mutable CurvilinearTrajectoryError  theCurvilinearError;
00153  
00154 };
00155 
00156 std::ostream& operator<<(std::ostream& os, const FreeTrajectoryState& fts);
00157 
00158 #endif