CMS 3D CMS Logo

CMSSW_4_4_3_patch1/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     theCartesianErrorValid(false),
00038     theCurvilinearErrorValid(false) {};
00039 
00040   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters) :
00041     theGlobalParameters(aGlobalParameters),
00042     theCartesianErrorValid(false),
00043     theCurvilinearErrorValid(false)
00044   {
00045   }
00046 
00047   FreeTrajectoryState(const GlobalPoint& aX,
00048                       const GlobalVector& aP,
00049                       TrackCharge aCharge, 
00050                       const MagneticField* fieldProvider) :
00051     theGlobalParameters(aX, aP, aCharge, fieldProvider),
00052     theCartesianErrorValid(false),
00053     theCurvilinearErrorValid(false)
00054   {
00055   }
00056 
00057   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00058                       const CartesianTrajectoryError& aCartesianError) :
00059     theGlobalParameters(aGlobalParameters),
00060     theCartesianError(aCartesianError),
00061     theCartesianErrorValid(true),
00062     theCurvilinearErrorValid(false)
00063   {
00064   }
00065   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00066                       const CurvilinearTrajectoryError& aCurvilinearError) :
00067     theGlobalParameters(aGlobalParameters),
00068     theCurvilinearError(aCurvilinearError),
00069     theCartesianErrorValid(false),
00070     theCurvilinearErrorValid(true)
00071   {
00072   }
00073   FreeTrajectoryState(const GlobalTrajectoryParameters& aGlobalParameters,
00074                       const CartesianTrajectoryError& aCartesianError,
00075                       const CurvilinearTrajectoryError& aCurvilinearError) :
00076     theGlobalParameters(aGlobalParameters),
00077     theCartesianError(aCartesianError),
00078     theCurvilinearError(aCurvilinearError),
00079     theCartesianErrorValid(true),
00080     theCurvilinearErrorValid(true)
00081   {
00082   }
00083 // access
00084 // propagate access to parameters
00085   GlobalPoint position() const {
00086     return theGlobalParameters.position();
00087   }
00088   GlobalVector momentum() const {
00089     return theGlobalParameters.momentum();
00090   }
00091   TrackCharge charge() const {
00092     return theGlobalParameters.charge();
00093   }
00094   double signedInverseMomentum() const {
00095     return theGlobalParameters.signedInverseMomentum();
00096   }
00097   double transverseCurvature() const {
00098     return theGlobalParameters.transverseCurvature();
00099   }
00100 
00101 // direct access
00102   bool hasCartesianError() const {return theCartesianErrorValid;}
00103   bool hasCurvilinearError() const {return theCurvilinearErrorValid;}
00104 
00105   bool hasError() const {
00106     return theCurvilinearErrorValid || theCartesianErrorValid;
00107   }
00108 
00109 
00110   const GlobalTrajectoryParameters& parameters() const {
00111     return theGlobalParameters;
00112   }
00113   const CartesianTrajectoryError& cartesianError() const {
00114     if unlikely(!hasError()) missingError();
00115     if  unlikely(!theCartesianErrorValid)
00116       createCartesianError();
00117     return theCartesianError;
00118   }
00119   const CurvilinearTrajectoryError& curvilinearError() const {
00120     if  unlikely(!hasError()) missingError();
00121     if  unlikely(!theCurvilinearErrorValid)
00122       createCurvilinearError();
00123     return theCurvilinearError;
00124   }
00125 
00126   void rescaleError(double factor);
00127 
00128   void setCartesianError(const CartesianTrajectoryError &err) {
00129         theCartesianError = err; theCartesianErrorValid = true;
00130   }
00131   void setCartesianError(const AlgebraicSymMatrix66 &err) {
00132         theCartesianError = CartesianTrajectoryError(err); theCartesianErrorValid = true;
00133   }
00134   void setCurvilinearError(const CurvilinearTrajectoryError &err) {
00135         theCurvilinearError = err; theCurvilinearErrorValid = true;
00136   }
00137   void setCurvilinearError(const AlgebraicSymMatrix55 &err) {
00138         theCurvilinearError = CurvilinearTrajectoryError(err); theCurvilinearErrorValid = true;
00139   }
00140 // properties
00141   bool canReach(double radius) const;
00142 private:
00143 
00144 
00145   static void missingError(); // dso_internal;
00146 
00147 // convert curvilinear errors to cartesian
00148   void createCartesianError() const; // dso_internal;
00149 // convert cartesian errors to curvilinear
00150   void createCurvilinearError() const; // dso_internal;
00151 
00152 private:
00153 
00154   GlobalTrajectoryParameters  theGlobalParameters;
00155   mutable CartesianTrajectoryError    theCartesianError;
00156   mutable CurvilinearTrajectoryError  theCurvilinearError;
00157   mutable bool                        theCartesianErrorValid;
00158   mutable bool                        theCurvilinearErrorValid;
00159 
00160 };
00161 
00162 std::ostream& operator<<(std::ostream& os, const FreeTrajectoryState& fts);
00163 
00164 #endif