Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:22

0001 #ifndef TrackReco_TrackBase_h
0002 #define TrackReco_TrackBase_h
0003 /** \class reco::TrackBase TrackBase.h DataFormats/TrackReco/interface/TrackBase.h
0004  *
0005  * Common base class to all track types, including Muon fits.
0006  * Internally, the following information is stored: <BR>
0007  *   <DT> A reference position on the track: (vx,vy,vz) </DT>
0008  *   <DT> Momentum at this given reference point on track: (px,py,pz) </DT>
0009  *   <DT> 5D curvilinear covariance matrix from the track fit </DT>
0010  *   <DT> Charge </DT>
0011  *   <DT> Chi-square and number of degrees of freedom </DT>
0012  *   <DT> Summary information of the hit pattern </DT>
0013  *
0014  * For tracks reconstructed in the CMS Tracker, the reference position is the point of
0015  * closest approach to the centre of CMS. For muons, this is not necessarily true.
0016  *
0017  * Parameters associated to the 5D curvilinear covariance matrix: <BR>
0018  * <B> (qoverp, lambda, phi, dxy, dsz) </B><BR>
0019  * defined as:  <BR>
0020  *   <DT> qoverp = q / abs(p) = signed inverse of momentum [1/GeV] </DT>
0021  *   <DT> lambda = pi/2 - polar angle at the given point </DT>
0022  *   <DT> phi = azimuth angle at the given point </DT>
0023  *   <DT> dxy = -vx*sin(phi) + vy*cos(phi) [cm] </DT>
0024  *   <DT> dsz = vz*cos(lambda) - (vx*cos(phi)+vy*sin(phi))*sin(lambda) [cm] </DT>
0025  *
0026  * Geometrically, dxy is the signed distance in the XY plane between the
0027  * the straight line passing through (vx,vy) with azimuthal angle phi and
0028  * the point (0,0).<BR>
0029  * The dsz parameter is the signed distance in the SZ plane between the
0030  * the straight line passing through (vx,vy,vz) with angles (phi, lambda) and
0031  * the point (s=0,z=0). The S axis is defined by the projection of the
0032  * straight line onto the XY plane. The convention is to assign the S
0033  * coordinate for (vx,vy) as the value vx*cos(phi)+vy*sin(phi). This value is
0034  * zero when (vx,vy) is the point of minimum transverse distance to (0,0).
0035  *
0036  * Note that dxy and dsz provide sensible estimates of the distance from
0037  * the true particle trajectory to (0,0,0) ONLY in two cases:<BR>
0038  *   <DT> When (vx,vy,vz) already correspond to the point of minimum transverse
0039  *   distance to (0,0,0) or is close to it (so that the differences
0040  *   between considering the exact trajectory or a straight line in this range
0041  *   are negligible). This is usually true for Tracker tracks. </DT>
0042  *   <DT> When the track has infinite or extremely high momentum </DT>
0043  *
0044  * More details about this parametrization are provided in the following document: <BR>
0045  * <a href="http://cms.cern.ch/iCMS/jsp/openfile.jsp?type=NOTE&year=2006&files=NOTE2006_001.pdf">A. Strandlie, W. Wittek, "Propagation of Covariance Matrices...", CMS Note 2006/001</a> <BR>
0046  *
0047  * \author Thomas Speer, Luca Lista, Pascal Vanlaer, Juan Alcaraz
0048  *
0049  */
0050 
0051 #include "DataFormats/TrackReco/interface/HitPattern.h"
0052 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0053 #include "DataFormats/Math/interface/Vector.h"
0054 #include "DataFormats/Math/interface/Error.h"
0055 #include "DataFormats/Math/interface/Vector3D.h"
0056 #include "DataFormats/Math/interface/Point3D.h"
0057 #include "DataFormats/Math/interface/Error.h"
0058 #include <bitset>
0059 
0060 namespace reco {
0061 
0062   class TrackBase {
0063   public:
0064     /// parameter dimension
0065     enum { dimension = 5 };
0066 
0067     /// error matrix size
0068     enum { covarianceSize = dimension * (dimension + 1) / 2 };
0069 
0070     /// parameter vector
0071     typedef math::Vector<dimension>::type ParameterVector;
0072 
0073     /// 5 parameter covariance matrix
0074     typedef math::Error<dimension>::type CovarianceMatrix;
0075 
0076     /// spatial vector
0077     typedef math::XYZVector Vector;
0078 
0079     /// point in the space
0080     typedef math::XYZPoint Point;
0081 
0082     /// enumerator provided indices to the five parameters
0083     enum { i_qoverp = 0, i_lambda, i_phi, i_dxy, i_dsz };
0084 
0085     /// index type
0086     typedef unsigned int index;
0087 
0088     /// track algorithm
0089     enum TrackAlgorithm {
0090       undefAlgorithm = 0,
0091       ctf = 1,
0092       duplicateMerge = 2,
0093       cosmics = 3,
0094       initialStep = 4,
0095       lowPtTripletStep = 5,
0096       pixelPairStep = 6,
0097       detachedTripletStep = 7,
0098       mixedTripletStep = 8,
0099       pixelLessStep = 9,
0100       tobTecStep = 10,
0101       jetCoreRegionalStep = 11,
0102       conversionStep = 12,
0103       muonSeededStepInOut = 13,
0104       muonSeededStepOutIn = 14,
0105       outInEcalSeededConv = 15,
0106       inOutEcalSeededConv = 16,
0107       nuclInter = 17,
0108       standAloneMuon = 18,
0109       globalMuon = 19,
0110       cosmicStandAloneMuon = 20,
0111       cosmicGlobalMuon = 21,
0112       // Phase1
0113       highPtTripletStep = 22,
0114       lowPtQuadStep = 23,
0115       detachedQuadStep = 24,
0116       displacedGeneralStep = 25,
0117       displacedRegionalStep = 26,
0118       bTagGhostTracks = 27,
0119       beamhalo = 28,
0120       gsf = 29,
0121       // HLT algo name
0122       hltPixel = 30,
0123       // steps used by PF
0124       hltIter0 = 31,
0125       hltIter1 = 32,
0126       hltIter2 = 33,
0127       hltIter3 = 34,
0128       hltIter4 = 35,
0129       // steps used by all other objects @HLT
0130       hltIterX = 36,
0131       // steps used by HI muon regional iterative tracking
0132       hiRegitMuInitialStep = 37,
0133       hiRegitMuLowPtTripletStep = 38,
0134       hiRegitMuPixelPairStep = 39,
0135       hiRegitMuDetachedTripletStep = 40,
0136       hiRegitMuMixedTripletStep = 41,
0137       hiRegitMuPixelLessStep = 42,
0138       hiRegitMuTobTecStep = 43,
0139       hiRegitMuMuonSeededStepInOut = 44,
0140       hiRegitMuMuonSeededStepOutIn = 45,
0141       algoSize = 46
0142     };
0143 
0144     /// algo mask
0145     typedef std::bitset<algoSize> AlgoMask;
0146 
0147     static const std::string algoNames[];
0148 
0149     /// track quality
0150     enum TrackQuality {
0151       undefQuality = -1,
0152       loose = 0,
0153       tight = 1,
0154       highPurity = 2,
0155       confirmed = 3,      // means found by more than one iteration
0156       goodIterative = 4,  // meaningless
0157       looseSetWithPV = 5,
0158       highPuritySetWithPV = 6,
0159       discarded = 7,  // because a better track found. kept in the collection for reference....
0160       qualitySize = 8
0161     };
0162 
0163     static const std::string qualityNames[];
0164 
0165     /// default constructor
0166     TrackBase();
0167 
0168     /// constructor from fit parameters and error matrix
0169     TrackBase(double chi2,
0170               double ndof,
0171               const Point &vertex,
0172               const Vector &momentum,
0173               int charge,
0174               const CovarianceMatrix &cov,
0175               TrackAlgorithm = undefAlgorithm,
0176               TrackQuality quality = undefQuality,
0177               signed char nloops = 0,
0178               uint8_t stopReason = 0,
0179               float t0 = 0.f,
0180               float beta = 0.f,
0181               float covt0t0 = -1.f,
0182               float covbetabeta = -1.f);
0183 
0184     /// virtual destructor
0185     virtual ~TrackBase();
0186 
0187     /// return true if timing measurement is usable
0188     bool isTimeOk() const { return covt0t0_ > 0.f; }
0189 
0190     /// chi-squared of the fit
0191     double chi2() const;
0192 
0193     /// number of degrees of freedom of the fit
0194     double ndof() const;
0195 
0196     /// chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
0197     double normalizedChi2() const;
0198 
0199     /// track electric charge
0200     int charge() const;
0201 
0202     /// q / p
0203     double qoverp() const;
0204 
0205     /// polar angle
0206     double theta() const;
0207 
0208     /// Lambda angle
0209     double lambda() const;
0210 
0211     /// dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close to (0,0,0): see parametrization definition above for details). See also function dxy(myBeamSpot).
0212     double dxy() const;
0213 
0214     /// dxy parameter in perigee convention (d0 = -dxy)
0215     double d0() const;
0216 
0217     /// dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from  (0,0,0): see parametrization definition above for details)
0218     double dsz() const;
0219 
0220     /// dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to (0,0,0). See also function dz(myBeamSpot)
0221     double dz() const;
0222 
0223     /// momentum vector magnitude square
0224     double p2() const;
0225 
0226     /// momentum vector magnitude
0227     double p() const;
0228 
0229     /// track transverse momentum square
0230     double pt2() const;
0231 
0232     /// track transverse momentum
0233     double pt() const;
0234 
0235     /// x coordinate of momentum vector
0236     double px() const;
0237 
0238     /// y coordinate of momentum vector
0239     double py() const;
0240 
0241     /// z coordinate of momentum vector
0242     double pz() const;
0243 
0244     /// azimuthal angle of momentum vector
0245     double phi() const;
0246 
0247     /// pseudorapidity of momentum vector
0248     double eta() const;
0249 
0250     /// x coordinate of the reference point on track
0251     double vx() const;
0252 
0253     /// y coordinate of the reference point on track
0254     double vy() const;
0255 
0256     /// z coordinate of the reference point on track
0257     double vz() const;
0258 
0259     /// track momentum vector
0260     const Vector &momentum() const;
0261 
0262     /// Reference point on the track
0263     const Point &referencePoint() const;
0264 
0265     /// time at the reference point
0266     double t0() const;
0267 
0268     /// velocity at the reference point in natural units
0269     double beta() const;
0270 
0271     /// reference point on the track. This method is DEPRECATED, please use referencePoint() instead
0272     const Point &vertex() const;
0273     //__attribute__((deprecated("This method is DEPRECATED, please use referencePoint() instead.")));
0274 
0275     /// dxy parameter with respect to a user-given beamSpot  (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.
0276     double dxy(const Point &myBeamSpot) const;
0277 
0278     /// dxy parameter with respect to the beamSpot taking into account the beamspot slopes (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.
0279     double dxy(const BeamSpot &theBeamSpot) const;
0280 
0281     /// dsz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the distance in the S-Z plane to the beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.
0282     double dsz(const Point &myBeamSpot) const;
0283 
0284     /// dz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the track z0, if the beamSpot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.
0285     double dz(const Point &myBeamSpot) const;
0286 
0287     /// Track parameters with one-to-one correspondence to the covariance matrix
0288     ParameterVector parameters() const;
0289 
0290     /// return track covariance matrix
0291     CovarianceMatrix covariance() const;
0292 
0293     /// i-th parameter ( i = 0, ... 4 )
0294     double parameter(int i) const;
0295 
0296     /// (i,j)-th element of covariance matrix (i, j = 0, ... 4)
0297     double covariance(int i, int j) const;
0298 
0299     /// error on t0
0300     double covt0t0() const;
0301 
0302     /// error on beta
0303     double covBetaBeta() const;
0304 
0305     /// error on specified element
0306     double error(int i) const;
0307 
0308     /// error on signed transverse curvature
0309     double qoverpError() const;
0310 
0311     /// error on Pt (set to 1000**2 TeV**2 if charge==0 for safety)
0312     double ptError2() const;
0313 
0314     /// error on Pt (set to 1000 TeV if charge==0 for safety)
0315     double ptError() const;
0316 
0317     /// error on theta
0318     double thetaError() const;
0319 
0320     /// error on lambda
0321     double lambdaError() const;
0322 
0323     /// error on eta
0324     double etaError() const;
0325 
0326     /// error on phi
0327     double phiError() const;
0328 
0329     /// error on dxy
0330     double dxyError() const;
0331 
0332     /// error on d0
0333     double d0Error() const;
0334 
0335     /// error on dsz
0336     double dszError() const;
0337 
0338     /// error on dz
0339     double dzError() const;
0340 
0341     /// error on t0
0342     double t0Error() const;
0343 
0344     /// error on beta
0345     double betaError() const;
0346 
0347     /// error on dxy with respect to a user-given reference point + uncertainty (i.e. reco::Vertex position)
0348     double dxyError(Point const &vtx, math::Error<3>::type const &vertexCov) const;
0349 
0350     /// error on dxy with respect to a user-given beamspot
0351     double dxyError(const BeamSpot &theBeamSpot) const;
0352 
0353     /// fill SMatrix
0354     CovarianceMatrix &fill(CovarianceMatrix &v) const;
0355 
0356     /// covariance matrix index in array
0357     static index covIndex(index i, index j);
0358 
0359     /// Access the hit pattern, indicating in which Tracker layers the track has hits.
0360     const HitPattern &hitPattern() const;
0361 
0362     /// number of valid hits found
0363     unsigned short numberOfValidHits() const;
0364 
0365     /// number of cases where track crossed a layer without getting a hit.
0366     unsigned short numberOfLostHits() const;
0367 
0368     /// number of hits expected from inner track extrapolation but missing
0369     int missingInnerHits() const;
0370 
0371     /// number of hits expected from outer track extrapolation but missing
0372     int missingOuterHits() const;
0373 
0374     /// fraction of valid hits on the track
0375     double validFraction() const;
0376 
0377     /// append hit patterns from vector of hit references
0378     template <typename C>
0379     bool appendHits(const C &c, const TrackerTopology &ttopo);
0380 
0381     template <typename I>
0382     bool appendHits(const I &begin, const I &end, const TrackerTopology &ttopo);
0383 
0384     /// append a single hit to the HitPattern
0385     bool appendHitPattern(const TrackingRecHit &hit, const TrackerTopology &ttopo);
0386     bool appendHitPattern(const DetId &id, TrackingRecHit::Type hitType, const TrackerTopology &ttopo);
0387 
0388     /**
0389      * These are meant to be used only in cases where the an
0390      * already-packed hit information is re-interpreted in terms of
0391      * HitPattern (i.e. MiniAOD PackedCandidate, and the IO rule for
0392      * reading old versions of HitPattern)
0393      */
0394     bool appendTrackerHitPattern(uint16_t subdet, uint16_t layer, uint16_t stereo, TrackingRecHit::Type hitType);
0395     bool appendHitPattern(const uint16_t pattern, TrackingRecHit::Type hitType);
0396 
0397     /**
0398      * This is meant to be used only in cases where the an
0399      * already-packed hit information is re-interpreted in terms of
0400      * HitPattern (i.e. the IO rule for reading old versions of
0401      * HitPattern)
0402      */
0403     bool appendMuonHitPattern(const DetId &id, TrackingRecHit::Type hitType);
0404 
0405     /// Sets HitPattern as empty
0406     void resetHitPattern();
0407 
0408     ///Track algorithm
0409     void setAlgorithm(const TrackAlgorithm a);
0410 
0411     void setOriginalAlgorithm(const TrackAlgorithm a);
0412 
0413     void setAlgoMask(AlgoMask a) { algoMask_ = a; }
0414 
0415     AlgoMask algoMask() const { return algoMask_; }
0416     unsigned long long algoMaskUL() const { return algoMask().to_ullong(); }
0417     bool isAlgoInMask(TrackAlgorithm a) const { return algoMask()[a]; }
0418 
0419     TrackAlgorithm algo() const;
0420     TrackAlgorithm originalAlgo() const;
0421 
0422     std::string algoName() const;
0423 
0424     static std::string algoName(TrackAlgorithm);
0425 
0426     static TrackAlgorithm algoByName(const std::string &name);
0427 
0428     ///Track quality
0429     bool quality(const TrackQuality) const;
0430 
0431     void setQuality(const TrackQuality);
0432 
0433     static std::string qualityName(TrackQuality);
0434 
0435     static TrackQuality qualityByName(const std::string &name);
0436 
0437     int qualityMask() const;
0438 
0439     void setQualityMask(int qualMask);
0440 
0441     void setNLoops(signed char value);
0442 
0443     bool isLooper() const;
0444 
0445     signed char nLoops() const;
0446 
0447     void setStopReason(uint8_t value) { stopReason_ = value; }
0448 
0449     uint8_t stopReason() const { return stopReason_; }
0450 
0451   private:
0452     /// hit pattern
0453     HitPattern hitPattern_;
0454 
0455     /// perigee 5x5 covariance matrix
0456     float covariance_[covarianceSize];
0457 
0458     /// errors for time and velocity (separate from cov for now)
0459     float covt0t0_, covbetabeta_;
0460 
0461     /// chi-squared
0462     float chi2_;
0463 
0464     /// innermost (reference) point on track
0465     Point vertex_;
0466 
0467     /// time at the reference point on track
0468     float t0_;
0469 
0470     /// momentum vector at innermost point
0471     Vector momentum_;
0472 
0473     /// norm of the particle velocity at innermost point on track
0474     /// can multiply by momentum_.Unit() to get velocity vector
0475     float beta_;
0476 
0477     /// algo mask, bit set for the algo where it was reconstructed + each algo a track was found overlapping by the listmerger
0478     std::bitset<algoSize> algoMask_;
0479 
0480     /// number of degrees of freedom
0481     float ndof_;
0482 
0483     /// electric charge
0484     char charge_;
0485 
0486     /// track algorithm
0487     uint8_t algorithm_;
0488 
0489     /// track algorithm
0490     uint8_t originalAlgorithm_;
0491 
0492     /// track quality
0493     uint8_t quality_;
0494 
0495     /// number of loops made during the building of the trajectory of a looper particle
0496     // I use signed char because I don't expect more than 128 loops and I could use a negative value for a special purpose.
0497     signed char nLoops_;
0498 
0499     /// Stop Reason
0500     uint8_t stopReason_;
0501   };
0502 
0503   //  Access the hit pattern, indicating in which Tracker layers the track has hits.
0504   inline const HitPattern &TrackBase::hitPattern() const { return hitPattern_; }
0505 
0506   inline bool TrackBase::appendHitPattern(const DetId &id, TrackingRecHit::Type hitType, const TrackerTopology &ttopo) {
0507     return hitPattern_.appendHit(id, hitType, ttopo);
0508   }
0509 
0510   inline bool TrackBase::appendHitPattern(const TrackingRecHit &hit, const TrackerTopology &ttopo) {
0511     return hitPattern_.appendHit(hit, ttopo);
0512   }
0513 
0514   inline bool TrackBase::appendTrackerHitPattern(uint16_t subdet,
0515                                                  uint16_t layer,
0516                                                  uint16_t stereo,
0517                                                  TrackingRecHit::Type hitType) {
0518     return hitPattern_.appendTrackerHit(subdet, layer, stereo, hitType);
0519   }
0520 
0521   inline bool TrackBase::appendHitPattern(uint16_t pattern, TrackingRecHit::Type hitType) {
0522     return hitPattern_.appendHit(pattern, hitType);
0523   }
0524 
0525   inline bool TrackBase::appendMuonHitPattern(const DetId &id, TrackingRecHit::Type hitType) {
0526     return hitPattern_.appendMuonHit(id, hitType);
0527   }
0528 
0529   inline void TrackBase::resetHitPattern() { hitPattern_.clear(); }
0530 
0531   template <typename I>
0532   bool TrackBase::appendHits(const I &begin, const I &end, const TrackerTopology &ttopo) {
0533     return hitPattern_.appendHits(begin, end, ttopo);
0534   }
0535 
0536   template <typename C>
0537   bool TrackBase::appendHits(const C &c, const TrackerTopology &ttopo) {
0538     return hitPattern_.appendHits(c.begin(), c.end(), ttopo);
0539   }
0540 
0541   inline TrackBase::index TrackBase::covIndex(index i, index j) {
0542     int a = (i <= j ? i : j);
0543     int b = (i <= j ? j : i);
0544     return b * (b + 1) / 2 + a;
0545   }
0546 
0547   inline TrackBase::TrackAlgorithm TrackBase::algo() const { return (TrackAlgorithm)(algorithm_); }
0548   inline TrackBase::TrackAlgorithm TrackBase::originalAlgo() const { return (TrackAlgorithm)(originalAlgorithm_); }
0549 
0550   inline std::string TrackBase::algoName() const { return TrackBase::algoName(algo()); }
0551 
0552   inline bool TrackBase::quality(const TrackBase::TrackQuality q) const {
0553     switch (q) {
0554       case undefQuality:
0555         return quality_ == 0;
0556       case goodIterative:
0557         return (quality_ & (1 << TrackBase::highPurity)) >> TrackBase::highPurity;
0558       default:
0559         return (quality_ & (1 << q)) >> q;
0560     }
0561     return false;
0562   }
0563 
0564   inline void TrackBase::setQuality(const TrackBase::TrackQuality q) {
0565     if (q == undefQuality) {
0566       quality_ = 0;
0567     } else {
0568       quality_ |= (1 << q);
0569     }
0570   }
0571 
0572   inline std::string TrackBase::qualityName(TrackQuality q) {
0573     if (int(q) < int(qualitySize) && int(q) >= 0) {
0574       return qualityNames[int(q)];
0575     }
0576     return "undefQuality";
0577   }
0578 
0579   inline std::string TrackBase::algoName(TrackAlgorithm a) {
0580     if (int(a) < int(algoSize) && int(a) > 0) {
0581       return algoNames[int(a)];
0582     }
0583     return "undefAlgorithm";
0584   }
0585 
0586   // chi-squared of the fit
0587   inline double TrackBase::chi2() const { return chi2_; }
0588 
0589   // number of degrees of freedom of the fit
0590   inline double TrackBase::ndof() const { return ndof_; }
0591 
0592   // chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
0593   inline double TrackBase::normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
0594 
0595   // track electric charge
0596   inline int TrackBase::charge() const { return charge_; }
0597 
0598   // q / p
0599   inline double TrackBase::qoverp() const { return charge() / p(); }
0600 
0601   // polar angle
0602   inline double TrackBase::theta() const { return momentum_.theta(); }
0603 
0604   // Lambda angle
0605   inline double TrackBase::lambda() const { return M_PI_2 - momentum_.theta(); }
0606 
0607   // dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close to (0,0,0): see parametrization definition above for details). See also function dxy(myBeamSpot) below.
0608   inline double TrackBase::dxy() const { return (-vx() * py() + vy() * px()) / pt(); }
0609 
0610   // dxy parameter in perigee convention (d0 = -dxy)
0611   inline double TrackBase::d0() const { return -dxy(); }
0612 
0613   // dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0,0,0): see parametrization definition above for details)
0614   inline double TrackBase::dsz() const {
0615     const auto thept = pt();
0616     const auto thepinv = 1 / p();
0617     const auto theptoverp = thept * thepinv;
0618     return vz() * theptoverp - (vx() * px() + vy() * py()) / thept * pz() * thepinv;
0619   }
0620 
0621   // dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to (0,0,0). See also function dz(myBeamSpot) below.
0622   inline double TrackBase::dz() const {
0623     const auto thept2inv = 1 / pt2();
0624     return vz() - (vx() * px() + vy() * py()) * pz() * thept2inv;
0625   }
0626 
0627   // momentum vector magnitude square
0628   inline double TrackBase::p2() const { return momentum_.Mag2(); }
0629 
0630   // momentum vector magnitude
0631   inline double TrackBase::p() const { return sqrt(p2()); }
0632 
0633   // track transverse momentum square
0634   inline double TrackBase::pt2() const { return momentum_.Perp2(); }
0635 
0636   // track transverse momentum
0637   inline double TrackBase::pt() const { return sqrt(pt2()); }
0638 
0639   // x coordinate of momentum vector
0640   inline double TrackBase::px() const { return momentum_.x(); }
0641 
0642   // y coordinate of momentum vector
0643   inline double TrackBase::py() const { return momentum_.y(); }
0644 
0645   // z coordinate of momentum vector
0646   inline double TrackBase::pz() const { return momentum_.z(); }
0647 
0648   // azimuthal angle of momentum vector
0649   inline double TrackBase::phi() const { return momentum_.Phi(); }
0650 
0651   // pseudorapidity of momentum vector
0652   inline double TrackBase::eta() const { return momentum_.Eta(); }
0653 
0654   // x coordinate of the reference point on track
0655   inline double TrackBase::vx() const { return vertex_.x(); }
0656 
0657   // y coordinate of the reference point on track
0658   inline double TrackBase::vy() const { return vertex_.y(); }
0659 
0660   // z coordinate of the reference point on track
0661   inline double TrackBase::vz() const { return vertex_.z(); }
0662 
0663   // track momentum vector
0664   inline const TrackBase::Vector &TrackBase::momentum() const { return momentum_; }
0665 
0666   // Reference point on the track
0667   inline const TrackBase::Point &TrackBase::referencePoint() const { return vertex_; }
0668 
0669   // Time at the reference point on the track
0670   inline double TrackBase::t0() const { return t0_; }
0671 
0672   // Velocity at the reference point on the track in natural units
0673   inline double TrackBase::beta() const { return beta_; }
0674 
0675   // reference point on the track. This method is DEPRECATED, please use referencePoint() instead
0676   inline const TrackBase::Point &TrackBase::vertex() const { return vertex_; }
0677 
0678   // dxy parameter with respect to a user-given beamSpot
0679   // (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved).
0680   // This is a good approximation for Tracker tracks.
0681   inline double TrackBase::dxy(const Point &myBeamSpot) const {
0682     return (-(vx() - myBeamSpot.x()) * py() + (vy() - myBeamSpot.y()) * px()) / pt();
0683   }
0684 
0685   // dxy parameter with respect to the beamSpot taking into account the beamspot slopes
0686   // (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved).
0687   // This is a good approximation for Tracker tracks.
0688   inline double TrackBase::dxy(const BeamSpot &theBeamSpot) const { return dxy(theBeamSpot.position(vz())); }
0689 
0690   // dsz parameter with respect to a user-given beamSpot
0691   // (WARNING: this quantity can only be interpreted as the distance in the S-Z plane to the beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved).
0692   // This is a good approximation for Tracker tracks.
0693   inline double TrackBase::dsz(const Point &myBeamSpot) const {
0694     const auto thept = pt();
0695     const auto thepinv = 1 / p();
0696     const auto theptoverp = thept * thepinv;
0697     return (vz() - myBeamSpot.z()) * theptoverp -
0698            ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / thept * pz() * thepinv;
0699   }
0700 
0701   // dz parameter with respect to a user-given beamSpot
0702   // (WARNING: this quantity can only be interpreted as the track z0, if the beamSpot is reasonably close to the refPoint, since linear approximations are involved).
0703   // This is a good approximation for Tracker tracks.
0704   inline double TrackBase::dz(const Point &myBeamSpot) const {
0705     const auto theptinv2 = 1 / pt2();
0706     return (vz() - myBeamSpot.z()) -
0707            ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) * pz() * theptinv2;
0708   }
0709 
0710   // Track parameters with one-to-one correspondence to the covariance matrix
0711   inline TrackBase::ParameterVector TrackBase::parameters() const {
0712     return TrackBase::ParameterVector(qoverp(), lambda(), phi(), dxy(), dsz());
0713   }
0714 
0715   // return track covariance matrix
0716   inline TrackBase::CovarianceMatrix TrackBase::covariance() const {
0717     CovarianceMatrix m;
0718     fill(m);
0719     return m;
0720   }
0721 
0722   // i-th parameter ( i = 0, ... 4 )
0723   inline double TrackBase::parameter(int i) const { return parameters()[i]; }
0724 
0725   // (i,j)-th element of covariance matrix (i, j = 0, ... 4)
0726   inline double TrackBase::covariance(int i, int j) const { return covariance_[covIndex(i, j)]; }
0727 
0728   // error on specified element
0729   inline double TrackBase::error(int i) const { return sqrt(covariance_[covIndex(i, i)]); }
0730 
0731   // error on signed transverse curvature
0732   inline double TrackBase::qoverpError() const { return error(i_qoverp); }
0733 
0734   // error on Pt (set to 1000**2 TeV**2 if charge==0 for safety)
0735   inline double TrackBase::ptError2() const {
0736     const auto thecharge = charge();
0737 
0738     if (thecharge != 0) {
0739       const auto thept2 = pt2();
0740       const auto thep2 = p2();
0741       const auto thepz = pz();
0742       const auto ptimespt = sqrt(thep2 * thept2);
0743       const auto oneovercharge = 1 / thecharge;
0744 
0745       return thept2 * thep2 * oneovercharge * oneovercharge * covariance(i_qoverp, i_qoverp) +
0746              2 * ptimespt * oneovercharge * thepz * covariance(i_qoverp, i_lambda) +
0747              thepz * thepz * covariance(i_lambda, i_lambda);
0748     }
0749 
0750     return 1.e12;
0751   }
0752 
0753   // error on Pt (set to 1000 TeV if charge==0 for safety)
0754   inline double TrackBase::ptError() const { return sqrt(ptError2()); }
0755 
0756   // error on theta
0757   inline double TrackBase::thetaError() const { return error(i_lambda); }
0758 
0759   // error on lambda
0760   inline double TrackBase::lambdaError() const { return error(i_lambda); }
0761 
0762   // error on eta
0763   inline double TrackBase::etaError() const { return error(i_lambda) * sqrt(p2() / pt2()); }
0764 
0765   // error on phi
0766   inline double TrackBase::phiError() const { return error(i_phi); }
0767 
0768   // error on dxy
0769   inline double TrackBase::dxyError() const { return error(i_dxy); }
0770 
0771   // error on d0
0772   inline double TrackBase::d0Error() const { return error(i_dxy); }
0773 
0774   // error on dsz
0775   inline double TrackBase::dszError() const { return error(i_dsz); }
0776 
0777   // error on dz
0778   inline double TrackBase::dzError() const { return error(i_dsz) * sqrt(p2() / pt2()); }
0779 
0780   // covariance of t0
0781   inline double TrackBase::covt0t0() const { return covt0t0_; }
0782 
0783   // covariance of beta
0784   inline double TrackBase::covBetaBeta() const { return covbetabeta_; }
0785 
0786   // error on t0
0787   inline double TrackBase::t0Error() const { return std::sqrt(covt0t0_); }
0788 
0789   // error on beta
0790   inline double TrackBase::betaError() const { return std::sqrt(covbetabeta_); }
0791 
0792   // error on dxy with respect to a given beamspot
0793   inline double TrackBase::dxyError(const BeamSpot &theBeamSpot) const {
0794     return dxyError(theBeamSpot.position(vz()), theBeamSpot.rotatedCovariance3D());
0795   }
0796 
0797   // number of valid hits found
0798   inline unsigned short TrackBase::numberOfValidHits() const { return hitPattern_.numberOfValidHits(); }
0799 
0800   // number of cases where track crossed a layer without getting a hit.
0801   inline unsigned short TrackBase::numberOfLostHits() const {
0802     return hitPattern_.numberOfLostHits(HitPattern::TRACK_HITS);
0803   }
0804 
0805   // number of hits expected from inner track extrapolation but missing
0806   inline int TrackBase::missingInnerHits() const {
0807     return hitPattern_.numberOfLostHits(HitPattern::MISSING_INNER_HITS);
0808   }
0809 
0810   // number of hits expected from outer track extrapolation but missing
0811   inline int TrackBase::missingOuterHits() const {
0812     return hitPattern_.numberOfLostHits(HitPattern::MISSING_OUTER_HITS);
0813   }
0814 
0815   // fraction of valid hits on the track
0816   inline double TrackBase::validFraction() const {
0817     int valid = hitPattern_.numberOfValidTrackerHits();
0818     int lost = hitPattern_.numberOfLostTrackerHits(HitPattern::TRACK_HITS);
0819     int lostIn = hitPattern_.numberOfLostTrackerHits(HitPattern::MISSING_INNER_HITS);
0820     int lostOut = hitPattern_.numberOfLostTrackerHits(HitPattern::MISSING_OUTER_HITS);
0821 
0822     const auto tot = valid + lost + lostIn + lostOut;
0823 
0824     if (tot == 0) {
0825       return -1;
0826     }
0827 
0828     return valid / (double)(tot);
0829   }
0830 
0831   //Track algorithm
0832   inline void TrackBase::setAlgorithm(const TrackBase::TrackAlgorithm a) {
0833     algorithm_ = a;
0834     algoMask_.reset();
0835     setOriginalAlgorithm(a);
0836   }
0837 
0838   inline void TrackBase::setOriginalAlgorithm(const TrackBase::TrackAlgorithm a) {
0839     originalAlgorithm_ = a;
0840     algoMask_.set(a);
0841   }
0842 
0843   inline int TrackBase::qualityMask() const { return quality_; }
0844 
0845   inline void TrackBase::setQualityMask(int qualMask) { quality_ = qualMask; }
0846 
0847   inline void TrackBase::setNLoops(signed char value) { nLoops_ = value; }
0848 
0849   inline bool TrackBase::isLooper() const { return (nLoops_ > 0); }
0850 
0851   inline signed char TrackBase::nLoops() const { return nLoops_; }
0852 
0853 }  // namespace reco
0854 
0855 #endif