CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackBase.h
Go to the documentation of this file.
1 #ifndef TrackReco_TrackBase_h
2 #define TrackReco_TrackBase_h
3 
58 
59 namespace reco
60 {
61 
62 class TrackBase
63 {
64 
65 public:
67  enum { dimension = 5 };
68 
70  enum { covarianceSize = dimension * (dimension + 1) / 2 };
71 
74 
77 
80 
83 
85  enum {
86  i_qoverp = 0,
91  };
92 
94  typedef unsigned int index;
95 
98  undefAlgorithm = 0, ctf = 1, rs = 2, cosmics = 3,
111  nuclInter = 17,
115  beamhalo = 28,
116  gsf = 29,
117  // HLT algo name
118  hltPixel = 30,
119  // steps used by PF
120  hltIter0 = 31,
121  hltIter1 = 32,
122  hltIter2 = 33,
123  hltIter3 = 34,
124  hltIter4 = 35,
125  // steps used by all other objects @HLT
126  hltIterX = 36,
127  algoSize = 37
128  };
129 
130  static const std::string algoNames[];
131 
135  loose = 0,
136  tight = 1,
143  };
144 
145  static const std::string qualityNames[];
146 
148  TrackBase();
149 
151  TrackBase(double chi2, double ndof, const Point &vertex,
152  const Vector &momentum, int charge, const CovarianceMatrix &cov,
154  signed char nloops = 0);
155 
157  virtual ~TrackBase();
158 
160  double chi2() const;
161 
163  double ndof() const;
164 
166  double normalizedChi2() const;
167 
169  int charge() const;
170 
172  double qoverp() const;
173 
175  double theta() const;
176 
178  double lambda() const;
179 
181  double dxy() const;
182 
184  double d0() const;
185 
187  double dsz() const;
188 
190  double dz() const;
191 
193  double p() const;
194 
196  double pt() const;
197 
199  double px() const;
200 
202  double py() const;
203 
205  double pz() const;
206 
208  double phi() const;
209 
211  double eta() const;
212 
214  double vx() const;
215 
217  double vy() const;
218 
220  double vz() const;
221 
223  const Vector &momentum() const;
224 
226  const Point &referencePoint() const;
227 
229  const Point &vertex() const ;
230  //__attribute__((deprecated("This method is DEPRECATED, please use referencePoint() instead.")));
231 
233  double dxy(const Point &myBeamSpot) const;
234 
236  double dxy(const BeamSpot &theBeamSpot) const;
237 
239  double dsz(const Point &myBeamSpot) const;
240 
242  double dz(const Point &myBeamSpot) const;
243 
245  ParameterVector parameters() const;
246 
248  CovarianceMatrix covariance() const;
249 
251  double parameter(int i) const;
252 
254  double covariance(int i, int j) const;
255 
257  double error(int i) const;
258 
260  double qoverpError() const;
261 
263  double ptError() const;
264 
266  double thetaError() const;
267 
269  double lambdaError() const;
270 
272  double etaError() const;
273 
275  double phiError() const;
276 
278  double dxyError() const;
279 
281  double d0Error() const;
282 
284  double dszError() const;
285 
287  double dzError() const;
288 
290  CovarianceMatrix &fill(CovarianceMatrix &v) const;
291 
293  static index covIndex(index i, index j);
294 
296  const HitPattern &hitPattern() const;
297 
299  unsigned short numberOfValidHits() const;
300 
302  unsigned short numberOfLostHits() const;
303 
305  double validFraction() const;
306 
308  template<typename C>
309  bool appendHits(const C &c);
310 
311  template<typename I>
312  bool appendHits(const I &begin, const I &end);
313 
315  bool appendHitPattern(const TrackingRecHit &hit);
316  bool appendHitPattern(const DetId &id, TrackingRecHit::Type hitType);
317 
319  void resetHitPattern();
320 
322  void setAlgorithm(const TrackAlgorithm a, bool set = true);
323 
324  TrackAlgorithm algo() const ;
325 
326  std::string algoName() const;
327 
329 
330  static TrackAlgorithm algoByName(const std::string &name);
331 
333  bool quality(const TrackQuality) const;
334 
335  void setQuality(const TrackQuality, bool set = true);
336 
338 
340 
341  int qualityMask() const;
342 
343  void setQualityMask(int qualMask);
344 
345  void setNLoops(signed char value);
346 
347  bool isLooper() const;
348 
349  signed char nLoops() const;
350 
351 private:
354 
357 
359  float chi2_;
360 
362  Point vertex_;
363 
365  Vector momentum_;
366 
368  float ndof_;
369 
371  char charge_;
372 
374  uint8_t algorithm_;
375 
377  uint8_t quality_;
378 
380  signed char nLoops_; // I use signed char because I don't expect more than 128 loops and I could use a negative value for a special purpose.
381 };
382 
383 // Access the hit pattern, indicating in which Tracker layers the track has hits.
384 inline const HitPattern & TrackBase::hitPattern() const
385 {
386  return hitPattern_;
387 }
388 
390 {
391  return hitPattern_.appendHit(id, hitType);
392 }
393 
395 {
396  return hitPattern_.appendHit(hit);
397 }
398 
400 {
401  hitPattern_.clear();
402 }
403 
404 template<typename I>
405 bool TrackBase::appendHits(const I &begin, const I &end)
406 {
407  return hitPattern_.appendHits(begin, end);
408 }
409 
410 template<typename C>
412 {
413  return setHitPattern(c.begin(), c.end());
414 }
415 
417 {
418  int a = (i <= j ? i : j);
419  int b = (i <= j ? j : i);
420  return b * (b + 1) / 2 + a;
421 }
422 
424 {
425  return (TrackAlgorithm) algorithm_;
426 }
427 
429 {
430  // I'd like to do:
431  // return TrackBase::algoName(algorithm_);
432  // but I cannot define a const static function. Why???
433  switch (algorithm_) {
434  case undefAlgorithm:
435  return "undefAlgorithm";
436  break;
437  case ctf:
438  return "ctf";
439  break;
440  case rs:
441  return "rs";
442  break;
443  case cosmics:
444  return "cosmics";
445  break;
446  case beamhalo:
447  return "beamhalo";
448  break;
449  case initialStep:
450  return "initialStep";
451  break;
452  case lowPtTripletStep:
453  return "lowPtTripletStep";
454  break;
455  case pixelPairStep:
456  return "pixelPairStep";
457  break;
458  case detachedTripletStep:
459  return "detachedTripletStep";
460  break;
461  case mixedTripletStep:
462  return "mixedTripletStep";
463  break;
464  case pixelLessStep:
465  return "pixelLessStep";
466  break;
467  case tobTecStep:
468  return "tobTecStep";
469  break;
470  case jetCoreRegionalStep:
471  return "jetCoreRegionalStep";
472  break;
473  case conversionStep:
474  return "conversionStep";
475  break;
476  case muonSeededStepInOut:
477  return "muonSeededStepInOut";
478  break;
479  case muonSeededStepOutIn:
480  return "muonSeededStepOutIn";
481  break;
482  case outInEcalSeededConv:
483  return "outInEcalSeededConv";
484  break;
485  case inOutEcalSeededConv:
486  return "inOutEcalSeededConv";
487  break;
488  case nuclInter:
489  return "nuclInter";
490  break;
491  case standAloneMuon:
492  return "standAloneMuon";
493  break;
494  case globalMuon:
495  return "globalMuon";
496  break;
498  return "cosmicStandAloneMuon";
499  break;
500  case cosmicGlobalMuon:
501  return "cosmicGlobalMuon";
502  break;
503  case iter1LargeD0:
504  return "iter1LargeD0";
505  break;
506  case iter2LargeD0:
507  return "iter2LargeD0";
508  break;
509  case iter3LargeD0:
510  return "iter3LargeD0";
511  break;
512  case iter4LargeD0:
513  return "iter4LargeD0";
514  break;
515  case iter5LargeD0:
516  return "iter5LargeD0";
517  break;
518  case bTagGhostTracks:
519  return "bTagGhostTracks";
520  break;
521  case gsf:
522  return "gsf";
523  break;
524  case hltPixel :
525  return "hltPixel";
526  break;
527  case hltIter0 :
528  return "hltIter0";
529  break;
530  case hltIter1 :
531  return "hltIter1";
532  break;
533  case hltIter2 :
534  return "hltIter2";
535  break;
536  case hltIter3 :
537  return "hltIter3";
538  break;
539  case hltIter4 :
540  return "hltIter4";
541  break;
542  case hltIterX :
543  return "hltIterX";
544  break;
545  default:
546  return "undefAlgorithm";
547  break;
548  }
549 }
550 
552 {
553  switch (q) {
554  case undefQuality:
555  return (quality_ == 0);
556  case goodIterative:
557  return (((quality_ & (1 << TrackBase::confirmed)) >> TrackBase::confirmed) ||
559  default:
560  return (quality_ & (1 << q)) >> q;
561  }
562  return false;
563 }
564 
565 inline void TrackBase::setQuality(const TrackBase::TrackQuality q, bool set)
566 {
567  if (q == undefQuality) {
568  quality_ = 0;
569  } else {
570  //regular OR if setting value to true
571  if (set) {
572  quality_ |= (1 << q);
573  } else {
574  // doing "half-XOR" if unsetting value
575  quality_ &= (~(1 << q));
576  }
577  }
578 }
579 
581 {
582  if (int(q) < int(qualitySize) && int(q) >= 0) {
583  return qualityNames[int(q)];
584  }
585  return "undefQuality";
586 }
587 
589 {
590  if (int(a) < int(algoSize) && int(a) > 0) {
591  return algoNames[int(a)];
592  }
593  return "undefAlgorithm";
594 }
595 
596 // chi-squared of the fit
597 inline double TrackBase::chi2() const
598 {
599  return chi2_;
600 }
601 
602 // number of degrees of freedom of the fit
603 inline double TrackBase::ndof() const
604 {
605  return ndof_;
606 }
607 
608 // chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
609 inline double TrackBase::normalizedChi2() const
610 {
611  return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6;
612 }
613 
614 // track electric charge
615 inline int TrackBase::charge() const
616 {
617  return charge_;
618 }
619 
620 // q / p
621 inline double TrackBase::qoverp() const
622 {
623  return charge() / p();
624 }
625 
626 // polar angle
627 inline double TrackBase::theta() const
628 {
629  return momentum_.theta();
630 }
631 
632 // Lambda angle
633 inline double TrackBase::lambda() const
634 {
635  return M_PI_2 - momentum_.theta();
636 }
637 
638 // 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.
639 inline double TrackBase::dxy() const
640 {
641  return (-vx() * py() + vy() * px()) / pt();
642 }
643 
644 // dxy parameter in perigee convention (d0 = -dxy)
645 inline double TrackBase::d0() const
646 {
647  return -dxy();
648 }
649 
650 // 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)
651 inline double TrackBase::dsz() const
652 {
653  return vz() * pt() / p() - (vx() * px() + vy() * py()) / pt() * pz() / p();
654 }
655 
656 // 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.
657 inline double TrackBase::dz() const
658 {
659  return vz() - (vx() * px() + vy() * py()) / pt() * (pz() / pt());
660 }
661 
662 // momentum vector magnitude
663 inline double TrackBase::p() const
664 {
665  return momentum_.R();
666 }
667 
668 // track transverse momentum
669 inline double TrackBase::pt() const
670 {
671  return sqrt(momentum_.Perp2());
672 }
673 
674 // x coordinate of momentum vector
675 inline double TrackBase::px() const
676 {
677  return momentum_.x();
678 }
679 
680 // y coordinate of momentum vector
681 inline double TrackBase::py() const
682 {
683  return momentum_.y();
684 }
685 
686 // z coordinate of momentum vector
687 inline double TrackBase::pz() const
688 {
689  return momentum_.z();
690 }
691 
692 // azimuthal angle of momentum vector
693 inline double TrackBase::phi() const
694 {
695  return momentum_.Phi();
696 }
697 
698 // pseudorapidity of momentum vector
699 inline double TrackBase::eta() const
700 {
701  return momentum_.Eta();
702 }
703 
704 // x coordinate of the reference point on track
705 inline double TrackBase::vx() const
706 {
707  return vertex_.x();
708 }
709 
710 // y coordinate of the reference point on track
711 inline double TrackBase::vy() const
712 {
713  return vertex_.y();
714 }
715 
716 // z coordinate of the reference point on track
717 inline double TrackBase::vz() const
718 {
719  return vertex_.z();
720 }
721 
722 // track momentum vector
724 {
725  return momentum_;
726 }
727 
728 // Reference point on the track
730 {
731  return vertex_;
732 }
733 
734 // reference point on the track. This method is DEPRECATED, please use referencePoint() instead
735 inline const TrackBase::Point & TrackBase::vertex() const
736 {
737  return vertex_;
738 }
739 
740 // dxy parameter with respect to a user-given beamSpot
741 // (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).
742 // This is a good approximation for Tracker tracks.
743 inline double TrackBase::dxy(const Point &myBeamSpot) const
744 {
745  return (-(vx() - myBeamSpot.x()) * py() + (vy() - myBeamSpot.y()) * px()) / pt();
746 }
747 
748 // dxy parameter with respect to the beamSpot taking into account the beamspot slopes
749 // (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).
750 // This is a good approximation for Tracker tracks.
751 inline double TrackBase::dxy(const BeamSpot &theBeamSpot) const
752 {
753  return dxy(theBeamSpot.position(vz()));
754 }
755 
756 // dsz parameter with respect to a user-given beamSpot
757 // (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).
758 // This is a good approximation for Tracker tracks.
759 inline double TrackBase::dsz(const Point &myBeamSpot) const
760 {
761  return (vz() - myBeamSpot.z()) * pt() / p() - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / p();
762 }
763 
764 // dz parameter with respect to a user-given beamSpot
765 // (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).
766 // This is a good approximation for Tracker tracks.
767 inline double TrackBase::dz(const Point &myBeamSpot) const
768 {
769  return (vz() - myBeamSpot.z()) - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / pt();
770 }
771 
772 // Track parameters with one-to-one correspondence to the covariance matrix
774 {
775  return TrackBase::ParameterVector(qoverp(), lambda(), phi(), dxy(), dsz());
776 }
777 
778 // return track covariance matrix
780 {
782  fill(m);
783  return m;
784 }
785 
786 // i-th parameter ( i = 0, ... 4 )
787 inline double TrackBase::parameter(int i) const
788 {
789  return parameters()[i];
790 }
791 
792 // (i,j)-th element of covariance matrix (i, j = 0, ... 4)
793 inline double TrackBase::covariance(int i, int j) const
794 {
795  return covariance_[covIndex(i, j)];
796 }
797 
798 // error on specified element
799 inline double TrackBase::error(int i) const
800 {
801  return sqrt(covariance_[covIndex(i, i)]);
802 }
803 
804 // error on signed transverse curvature
805 inline double TrackBase::qoverpError() const
806 {
807  return error(i_qoverp);
808 }
809 
810 // error on Pt (set to 1000 TeV if charge==0 for safety)
811 inline double TrackBase::ptError() const
812 {
813  return (charge() != 0) ? sqrt(
814  pt() * pt() * p() * p() / charge() / charge() * covariance(i_qoverp, i_qoverp)
815  + 2 * pt() * p() / charge() * pz() * covariance(i_qoverp, i_lambda)
816  + pz() * pz() * covariance(i_lambda, i_lambda)) : 1.e6;
817 }
818 
819 // error on theta
820 inline double TrackBase::thetaError() const
821 {
822  return error(i_lambda);
823 }
824 
825 // error on lambda
826 inline double TrackBase::lambdaError() const
827 {
828  return error(i_lambda);
829 }
830 
831 // error on eta
832 inline double TrackBase::etaError() const
833 {
834  return error(i_lambda) * p() / pt();
835 }
836 
837 // error on phi
838 inline double TrackBase::phiError() const
839 {
840  return error(i_phi);
841 }
842 
843 // error on dxy
844 inline double TrackBase::dxyError() const
845 {
846  return error(i_dxy);
847 }
848 
849 // error on d0
850 inline double TrackBase::d0Error() const
851 {
852  return error(i_dxy);
853 }
854 
855 // error on dsz
856 inline double TrackBase::dszError() const
857 {
858  return error(i_dsz);
859 }
860 
861 // error on dz
862 inline double TrackBase::dzError() const
863 {
864  return error(i_dsz) * p() / pt();
865 }
866 
867 // number of valid hits found
868 inline unsigned short TrackBase::numberOfValidHits() const
869 {
871 }
872 
873 // number of cases where track crossed a layer without getting a hit.
874 inline unsigned short TrackBase::numberOfLostHits() const
875 {
877 }
878 
879 // fraction of valid hits on the track
880 inline double TrackBase::validFraction() const
881 {
886 
887  if ((valid + lost + lostIn + lostOut) == 0) {
888  return -1;
889  }
890 
891  return valid / (double)(valid + lost + lostIn + lostOut);
892 }
893 
894 //Track algorithm
896 {
897  if (set) {
898  algorithm_ = a;
899  } else {
901  }
902 }
903 
904 inline int TrackBase::qualityMask() const
905 {
906  return quality_;
907 }
908 
909 inline void TrackBase::setQualityMask(int qualMask)
910 {
911  quality_ = qualMask;
912 }
913 
914 inline void TrackBase::setNLoops(signed char value)
915 {
916  nLoops_ = value;
917 }
918 
919 inline bool TrackBase::isLooper() const
920 {
921  return (nLoops_ > 0);
922 }
923 
924 inline signed char TrackBase::nLoops() const
925 {
926  return nLoops_;
927 }
928 
929 } // namespace reco
930 
931 #endif
932 
double qoverp() const
q / p
Definition: TrackBase.h:621
double p() const
momentum vector magnitude
Definition: TrackBase.h:663
float chi2_
chi-squared
Definition: TrackBase.h:359
void setQualityMask(int qualMask)
Definition: TrackBase.h:909
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:729
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:850
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:580
unsigned int index
index type
Definition: TrackBase.h:94
bool isLooper() const
Definition: TrackBase.h:919
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:416
uint8_t quality_
track quality
Definition: TrackBase.h:377
double validFraction() const
fraction of valid hits on the track
Definition: TrackBase.h:880
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:645
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:609
TrackBase()
default constructor
Definition: TrackBase.cc:59
TrackQuality
track quality
Definition: TrackBase.h:133
double theta() const
polar angle
Definition: TrackBase.h:627
double dxyError() const
error on dxy
Definition: TrackBase.h:844
char charge_
electric charge
Definition: TrackBase.h:371
int numberOfValidHits() const
Definition: HitPattern.h:734
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:362
double etaError() const
error on eta
Definition: TrackBase.h:832
#define M_PI_2
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:693
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:874
ErrorD< N >::type type
Definition: Error.h:39
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:73
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:675
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:723
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
Definition: TrackBase.h:651
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:380
TrackAlgorithm
track algorithm
Definition: TrackBase.h:97
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:809
TrackAlgorithm algo() const
Definition: TrackBase.h:423
void setNLoops(signed char value)
Definition: TrackBase.h:914
bool appendHits(const C &c)
append hit patterns from vector of hit references
Definition: TrackBase.h:411
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:735
double dszError() const
error on dsz
Definition: TrackBase.h:856
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:699
void setQuality(const TrackQuality, bool set=true)
Definition: TrackBase.h:565
static const std::string qualityNames[]
Definition: TrackBase.h:145
fixed size vector
Definition: Vector.h:31
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:597
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:603
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:779
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:669
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:811
double phiError() const
error on phi
Definition: TrackBase.h:838
int qualityMask() const
Definition: TrackBase.h:904
int j
Definition: DBlmapReader.cc:9
bool appendHitPattern(const TrackingRecHit &hit)
append a single hit to the HitPattern
Definition: TrackBase.h:394
double lambda() const
Lambda angle.
Definition: TrackBase.h:633
const std::complex< double > I
Definition: I.h:8
double error(int i) const
error on specified element
Definition: TrackBase.h:799
bool appendHits(const I &begin, const I &end)
Definition: HitPattern.h:446
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:868
void resetHitPattern()
Sets HitPattern as empty.
Definition: TrackBase.h:399
#define end
Definition: vmac.h:37
math::XYZPoint Point
point in the space
Definition: TrackBase.h:82
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:687
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:787
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:805
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:657
double dzError() const
error on dz
Definition: TrackBase.h:862
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:717
Definition: DetId.h:18
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:108
virtual ~TrackBase()
virtual destructor
Definition: TrackBase.cc:98
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:384
signed char nLoops() const
Definition: TrackBase.h:924
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::string algoName() const
Definition: TrackBase.h:428
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:28
double b
Definition: hdecay.h:120
float ndof_
number of degrees of freedom
Definition: TrackBase.h:368
static const std::string algoNames[]
Definition: TrackBase.h:130
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:551
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:804
int numberOfValidTrackerHits() const
Definition: HitPattern.h:739
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:773
double lambdaError() const
error on lambda
Definition: TrackBase.h:826
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:711
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:356
#define begin
Definition: vmac.h:30
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:353
void setAlgorithm(const TrackAlgorithm a, bool set=true)
Track algorithm.
Definition: TrackBase.h:895
double a
Definition: hdecay.h:121
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:103
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:120
bool appendHit(const TrackingRecHit &hit)
Definition: HitPattern.cc:168
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:365
int charge() const
track electric charge
Definition: TrackBase.h:615
const Point & position() const
position
Definition: BeamSpot.h:62
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:374
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:639
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:79
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:681
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:705
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:76
double thetaError() const
error on theta
Definition: TrackBase.h:820