File indexing completed on 2024-04-06 12:04:16
0001 #ifndef GsfTrackReco_GsfTrack_h
0002 #define GsfTrackReco_GsfTrack_h
0003
0004
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
0015 enum { dimensionMode = 3 };
0016
0017 enum { covarianceSizeMode = dimensionMode * (dimensionMode + 1) / 2 };
0018
0019 typedef math::Vector<dimensionMode>::type ParameterVectorMode;
0020
0021 typedef math::Error<dimensionMode>::type CovarianceMatrixMode;
0022
0023 GsfTrack();
0024
0025
0026
0027 GsfTrack(double chi2, double ndof, const Point&, const Vector&, int charge, const CovarianceMatrix&);
0028
0029 void setGsfExtra(const GsfTrackExtraRef& ref) { gsfExtra_ = ref; }
0030
0031 const GsfTrackExtraRef& gsfExtra() const { return gsfExtra_; }
0032
0033
0034 void setMode(int chargeMode, const Vector& momentumMode, const CovarianceMatrixMode& covarianceMode);
0035
0036
0037 int chargeMode() const { return chargeMode_; }
0038
0039 double qoverpMode() const { return chargeMode() / pMode(); }
0040
0041 double thetaMode() const { return momentumMode_.theta(); }
0042
0043 double lambdaMode() const { return M_PI / 2 - momentumMode_.theta(); }
0044
0045 double pMode() const { return momentumMode_.R(); }
0046
0047 double ptMode() const { return sqrt(momentumMode_.Perp2()); }
0048
0049 double pxMode() const { return momentumMode_.x(); }
0050
0051 double pyMode() const { return momentumMode_.y(); }
0052
0053 double pzMode() const { return momentumMode_.z(); }
0054
0055 double phiMode() const { return momentumMode_.Phi(); }
0056
0057 double etaMode() const { return momentumMode_.Eta(); }
0058
0059
0060 const Vector& momentumMode() const { return momentumMode_; }
0061
0062
0063 ParameterVectorMode parametersMode() const { return ParameterVectorMode(qoverpMode(), lambdaMode(), phiMode()); }
0064
0065 CovarianceMatrixMode covarianceMode() const {
0066 CovarianceMatrixMode m;
0067 fill(m);
0068 return m;
0069 }
0070
0071
0072 double parameterMode(int i) const { return parametersMode()[i]; }
0073
0074 double covarianceMode(int i, int j) const { return covarianceMode_[covIndex(i, j)]; }
0075
0076 double errorMode(int i) const { return sqrt(covarianceMode_[covIndex(i, i)]); }
0077
0078
0079 double qoverpModeError() const { return errorMode(i_qoverp); }
0080
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
0090 double thetaModeError() const { return errorMode(i_lambda); }
0091
0092 double lambdaModeError() const { return errorMode(i_lambda); }
0093
0094 double etaModeError() const { return errorMode(i_lambda) * pMode() / ptMode(); }
0095
0096 double phiModeError() const { return errorMode(i_phi); }
0097
0098 private:
0099
0100 CovarianceMatrixMode& fill CMS_THREAD_SAFE(CovarianceMatrixMode& v) const;
0101
0102 private:
0103
0104 GsfTrackExtraRef gsfExtra_;
0105
0106 char chargeMode_;
0107
0108 Vector momentumMode_;
0109
0110 float covarianceMode_[covarianceSizeMode];
0111 };
0112
0113 }
0114
0115 #endif