Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:16

0001 #ifndef GsfTrackReco_GsfTrack_h
0002 #define GsfTrackReco_GsfTrack_h
0003 /** Extension of reco::Track for GSF. It contains
0004  * one additional Ref to a GsfTrackExtra object.
0005  */
0006 #include "DataFormats/TrackReco/interface/Track.h"
0007 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtraFwd.h"
0008 #include "FWCore/Utilities/interface/thread_safety_macros.h"
0009 
0010 namespace reco {
0011 
0012   class GsfTrack : public Track {
0013   public:
0014     /// parameter dimension mode
0015     enum { dimensionMode = 3 };
0016     /// error matrix size mode
0017     enum { covarianceSizeMode = dimensionMode * (dimensionMode + 1) / 2 };
0018     /// parameter vector (momentum part) from mode
0019     typedef math::Vector<dimensionMode>::type ParameterVectorMode;
0020     /// 3 parameter covariance matrix (momentum part) from mode
0021     typedef math::Error<dimensionMode>::type CovarianceMatrixMode;
0022     /// default constructor
0023     GsfTrack();
0024     /// constructor from fit parameters and error matrix
0025     /// notice that the reference point must be
0026     /// the point of closest approch to the beamline.
0027     GsfTrack(double chi2, double ndof, const Point&, const Vector&, int charge, const CovarianceMatrix&);
0028     /// set reference to GSF "extra" object
0029     void setGsfExtra(const GsfTrackExtraRef& ref) { gsfExtra_ = ref; }
0030     /// reference to "extra" object
0031     const GsfTrackExtraRef& gsfExtra() const { return gsfExtra_; }
0032 
0033     /// set mode parameters
0034     void setMode(int chargeMode, const Vector& momentumMode, const CovarianceMatrixMode& covarianceMode);
0035 
0036     /// track electric charge from mode
0037     int chargeMode() const { return chargeMode_; }
0038     /// q/p  from mode
0039     double qoverpMode() const { return chargeMode() / pMode(); }
0040     /// polar angle   from mode
0041     double thetaMode() const { return momentumMode_.theta(); }
0042     /// Lambda angle from mode
0043     double lambdaMode() const { return M_PI / 2 - momentumMode_.theta(); }
0044     /// momentum vector magnitude from mode
0045     double pMode() const { return momentumMode_.R(); }
0046     /// track transverse momentum from mode
0047     double ptMode() const { return sqrt(momentumMode_.Perp2()); }
0048     /// x coordinate of momentum vector from mode
0049     double pxMode() const { return momentumMode_.x(); }
0050     /// y coordinate of momentum vector from mode
0051     double pyMode() const { return momentumMode_.y(); }
0052     /// z coordinate of momentum vector from mode
0053     double pzMode() const { return momentumMode_.z(); }
0054     /// azimuthal angle of momentum vector from mode
0055     double phiMode() const { return momentumMode_.Phi(); }
0056     /// pseudorapidity of momentum vector from mode
0057     double etaMode() const { return momentumMode_.Eta(); }
0058 
0059     /// track momentum vector from mode
0060     const Vector& momentumMode() const { return momentumMode_; }
0061 
0062     /// Track parameters with one-to-one correspondence to the covariance matrix from mode
0063     ParameterVectorMode parametersMode() const { return ParameterVectorMode(qoverpMode(), lambdaMode(), phiMode()); }
0064     /// return track covariance matrix from mode
0065     CovarianceMatrixMode covarianceMode() const {
0066       CovarianceMatrixMode m;
0067       fill(m);
0068       return m;
0069     }
0070 
0071     /// i-th parameter ( i = 0, ... 2 ) from mode
0072     double parameterMode(int i) const { return parametersMode()[i]; }
0073     /// (i,j)-th element of covarianve matrix ( i, j = 0, ... 2 ) from mode
0074     double covarianceMode(int i, int j) const { return covarianceMode_[covIndex(i, j)]; }
0075     /// error on specified element from mode
0076     double errorMode(int i) const { return sqrt(covarianceMode_[covIndex(i, i)]); }
0077 
0078     /// error on signed transverse curvature from mode
0079     double qoverpModeError() const { return errorMode(i_qoverp); }
0080     /// error on Pt (set to 1000 TeV if charge==0 for safety) from mode
0081     double ptModeError() const {
0082       return (chargeMode() != 0)
0083                  ? sqrt(ptMode() * ptMode() * pMode() * pMode() / chargeMode() / chargeMode() *
0084                             covarianceMode(i_qoverp, i_qoverp) +
0085                         2 * ptMode() * pMode() / chargeMode() * pzMode() * covarianceMode(i_qoverp, i_lambda) +
0086                         pzMode() * pzMode() * covarianceMode(i_lambda, i_lambda))
0087                  : 1.e6;
0088     }
0089     /// error on theta from mode
0090     double thetaModeError() const { return errorMode(i_lambda); }
0091     /// error on lambda from mode
0092     double lambdaModeError() const { return errorMode(i_lambda); }
0093     /// error on eta from mode
0094     double etaModeError() const { return errorMode(i_lambda) * pMode() / ptMode(); }
0095     /// error on phi from mode
0096     double phiModeError() const { return errorMode(i_phi); }
0097 
0098   private:
0099     /// fill 3x3 SMatrix
0100     CovarianceMatrixMode& fill CMS_THREAD_SAFE(CovarianceMatrixMode& v) const;
0101 
0102   private:
0103     /// reference to GSF "extra" extension
0104     GsfTrackExtraRef gsfExtra_;
0105     /// electric charge from mode
0106     char chargeMode_;
0107     /// momentum vector from mode
0108     Vector momentumMode_;
0109     /// 3x3 momentum part of covariance (in q/p, lambda, phi)
0110     float covarianceMode_[covarianceSizeMode];
0111   };
0112 
0113 }  // namespace reco
0114 
0115 #endif