CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/DataFormats/GsfTrackReco/interface/GsfTrack.h

Go to the documentation of this file.
00001 #ifndef GsfTrackReco_GsfTrack_h
00002 #define GsfTrackReco_GsfTrack_h
00003 
00006 #include "DataFormats/TrackReco/interface/Track.h"
00007 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtraFwd.h" 
00008 
00009 namespace reco {
00010 
00011   class GsfTrack : public Track {
00012   public:
00014     enum { dimensionMode = 3 };
00016     enum { covarianceSizeMode = dimensionMode * ( dimensionMode + 1 ) / 2 };
00018     typedef math::Vector<dimensionMode>::type ParameterVectorMode;
00020     typedef math::Error<dimensionMode>::type CovarianceMatrixMode;
00022     GsfTrack();
00026     GsfTrack( double chi2, double ndof, const Point &, const Vector &, int charge,
00027               const CovarianceMatrix & );
00029     void setGsfExtra( const GsfTrackExtraRef & ref ) { gsfExtra_ = ref; }
00031     const GsfTrackExtraRef & gsfExtra() const { return gsfExtra_; }
00032 
00034     void setMode (int chargeMode, const Vector& momentumMode,
00035                   const CovarianceMatrixMode& covarianceMode);
00036 
00038     int chargeMode() const { return chargeMode_; }
00040     double qoverpMode() const { return chargeMode() / pMode(); }
00042     double thetaMode() const { return momentumMode_.theta(); }
00044     double lambdaMode() const { return M_PI/2 - momentumMode_.theta(); }
00046     double pMode() const { return momentumMode_.R(); }
00048     double ptMode() const { return sqrt( momentumMode_.Perp2() ); }
00050     double pxMode() const { return momentumMode_.x(); }
00052     double pyMode() const { return momentumMode_.y(); }
00054     double pzMode() const { return momentumMode_.z(); }
00056     double phiMode() const { return momentumMode_.Phi(); }
00058     double etaMode() const { return momentumMode_.Eta(); }
00059 
00061     const Vector & momentumMode() const { return momentumMode_; }
00062 
00064     ParameterVectorMode parametersMode() const { 
00065       return ParameterVectorMode(qoverpMode(),lambdaMode(),phiMode());
00066     }
00068     CovarianceMatrixMode covarianceMode() const { CovarianceMatrixMode m; fill( m ); return m; }
00069 
00071     double parameterMode(int i) const { return parametersMode()[i]; }
00073     double covarianceMode( int i, int j ) const { return covarianceMode_[ covIndex( i, j ) ]; }
00075     double errorMode( int i ) const { return sqrt( covarianceMode_[ covIndex( i, i ) ] ); }
00076 
00078     double qoverpModeError() const { return errorMode( i_qoverp ); }
00080     double ptModeError() const { 
00081       return (chargeMode()!=0) ?  sqrt( 
00082             ptMode()*ptMode()*pMode()*pMode()/chargeMode()/chargeMode() * covarianceMode(i_qoverp,i_qoverp)
00083           + 2*ptMode()*pMode()/chargeMode()*pzMode() * covarianceMode(i_qoverp,i_lambda)
00084           + pzMode()*pzMode() * covarianceMode(i_lambda,i_lambda) ) : 1.e6;
00085     }
00087     double thetaModeError() const { return errorMode( i_lambda ); }
00089     double lambdaModeError() const { return errorMode( i_lambda ); }
00091     double etaModeError() const { return errorMode( i_lambda ) * pMode()/ptMode(); }
00093     double phiModeError() const { return errorMode( i_phi ); }
00094 
00095   private:
00097     CovarianceMatrixMode & fill( CovarianceMatrixMode & v ) const;
00098 
00099 
00100   private:
00102     GsfTrackExtraRef gsfExtra_;
00104     char chargeMode_;
00106     Vector momentumMode_;
00108     float covarianceMode_[ covarianceSizeMode ];
00109 
00110   };
00111 
00112 }
00113 
00114 #endif