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 #include <bitset>
59 
60 namespace reco
61 {
62 
63 class TrackBase
64 {
65 
66 public:
68  enum { dimension = 5 };
69 
71  enum { covarianceSize = dimension * (dimension + 1) / 2 };
72 
75 
78 
81 
84 
86  enum {
87  i_qoverp = 0,
92  };
93 
95  typedef unsigned int index;
96 
97 
100  undefAlgorithm = 0, ctf = 1, rs = 2, cosmics = 3,
113  nuclInter = 17,
117  beamhalo = 28,
118  gsf = 29,
119  // HLT algo name
120  hltPixel = 30,
121  // steps used by PF
122  hltIter0 = 31,
123  hltIter1 = 32,
124  hltIter2 = 33,
125  hltIter3 = 34,
126  hltIter4 = 35,
127  // steps used by all other objects @HLT
128  hltIterX = 36,
129  algoSize = 37
130  };
131 
133  typedef std::bitset<algoSize> AlgoMask;
134 
135 
136  static const std::string algoNames[];
137 
141  loose = 0,
142  tight = 1,
144  confirmed = 3, // means found by more than one iteration
145  goodIterative = 4, // meaningless
148  discarded = 7, // because a better track found. kept in the collection for reference....
150  };
151 
152  static const std::string qualityNames[];
153 
155  TrackBase();
156 
158  TrackBase(double chi2, double ndof, const Point &vertex,
159  const Vector &momentum, int charge, const CovarianceMatrix &cov,
161  signed char nloops = 0);
162 
164  virtual ~TrackBase();
165 
167  double chi2() const;
168 
170  double ndof() const;
171 
173  double normalizedChi2() const;
174 
176  int charge() const;
177 
179  double qoverp() const;
180 
182  double theta() const;
183 
185  double lambda() const;
186 
188  double dxy() const;
189 
191  double d0() const;
192 
194  double dsz() const;
195 
197  double dz() const;
198 
200  double p() const;
201 
203  double pt() const;
204 
206  double px() const;
207 
209  double py() const;
210 
212  double pz() const;
213 
215  double phi() const;
216 
218  double eta() const;
219 
221  double vx() const;
222 
224  double vy() const;
225 
227  double vz() const;
228 
230  const Vector &momentum() const;
231 
233  const Point &referencePoint() const;
234 
236  const Point &vertex() const ;
237  //__attribute__((deprecated("This method is DEPRECATED, please use referencePoint() instead.")));
238 
240  double dxy(const Point &myBeamSpot) const;
241 
243  double dxy(const BeamSpot &theBeamSpot) const;
244 
246  double dsz(const Point &myBeamSpot) const;
247 
249  double dz(const Point &myBeamSpot) const;
250 
252  ParameterVector parameters() const;
253 
255  CovarianceMatrix covariance() const;
256 
258  double parameter(int i) const;
259 
261  double covariance(int i, int j) const;
262 
264  double error(int i) const;
265 
267  double qoverpError() const;
268 
270  double ptError() const;
271 
273  double thetaError() const;
274 
276  double lambdaError() const;
277 
279  double etaError() const;
280 
282  double phiError() const;
283 
285  double dxyError() const;
286 
288  double d0Error() const;
289 
291  double dszError() const;
292 
294  double dzError() const;
295 
297  CovarianceMatrix &fill(CovarianceMatrix &v) const;
298 
300  static index covIndex(index i, index j);
301 
303  const HitPattern &hitPattern() const;
304 
306  unsigned short numberOfValidHits() const;
307 
309  unsigned short numberOfLostHits() const;
310 
312  double validFraction() const;
313 
315  template<typename C>
316  bool appendHits(const C &c);
317 
318  template<typename I>
319  bool appendHits(const I &begin, const I &end);
320 
322  bool appendHitPattern(const TrackingRecHit &hit);
323  bool appendHitPattern(const DetId &id, TrackingRecHit::Type hitType);
324 
326  void resetHitPattern();
327 
329  void setAlgorithm(const TrackAlgorithm a);
330 
332 
333  void setAlgoMask(AlgoMask a) { algoMask_ = a;}
334 
335  AlgoMask algoMask() const { return algoMask_;}
336 #if ( !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__) ) || defined(__ROOTCLING__)
337  unsigned long long algoMaskUL() const { return algoMask().to_ullong();}
338 #endif
339  bool isAlgoInMask(TrackAlgorithm a) const {return algoMask()[a];}
340 
341 
342  TrackAlgorithm algo() const ;
343  TrackAlgorithm originalAlgo() const ;
344 
345 
346  std::string algoName() const;
347 
349 
350  static TrackAlgorithm algoByName(const std::string &name);
351 
353  bool quality(const TrackQuality) const;
354 
355  void setQuality(const TrackQuality);
356 
358 
360 
361  int qualityMask() const;
362 
363  void setQualityMask(int qualMask);
364 
365  void setNLoops(signed char value);
366 
367  bool isLooper() const;
368 
369  signed char nLoops() const;
370 
371 private:
374 
377 
379  float chi2_;
380 
382  Point vertex_;
383 
385  Vector momentum_;
386 
388  std::bitset<algoSize> algoMask_;
389 
391  float ndof_;
392 
394  char charge_;
395 
397  uint8_t algorithm_;
398 
401 
402 
404  uint8_t quality_;
405 
407  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.
408 };
409 
410 // Access the hit pattern, indicating in which Tracker layers the track has hits.
411 inline const HitPattern & TrackBase::hitPattern() const
412 {
413  return hitPattern_;
414 }
415 
417 {
418  return hitPattern_.appendHit(id, hitType);
419 }
420 
422 {
423  return hitPattern_.appendHit(hit);
424 }
425 
427 {
428  hitPattern_.clear();
429 }
430 
431 template<typename I>
432 bool TrackBase::appendHits(const I &begin, const I &end)
433 {
434  return hitPattern_.appendHits(begin, end);
435 }
436 
437 template<typename C>
439 {
440  return setHitPattern(c.begin(), c.end());
441 }
442 
444 {
445  int a = (i <= j ? i : j);
446  int b = (i <= j ? j : i);
447  return b * (b + 1) / 2 + a;
448 }
449 
451 {
452  return (TrackAlgorithm) (algorithm_);
453 }
455 {
457 }
458 
459 
460 
462 
464 {
465  switch (q) {
466  case undefQuality:
467  return quality_ == 0;
468  case goodIterative:
470  default:
471  return (quality_ & (1 << q)) >> q;
472  }
473  return false;
474 }
475 
477 {
478  if (q == undefQuality) {
479  quality_ = 0;
480  } else {
481  quality_ |= (1 << q);
482  }
483 }
484 
486 {
487  if (int(q) < int(qualitySize) && int(q) >= 0) {
488  return qualityNames[int(q)];
489  }
490  return "undefQuality";
491 }
492 
494 {
495  if (int(a) < int(algoSize) && int(a) > 0) {
496  return algoNames[int(a)];
497  }
498  return "undefAlgorithm";
499 }
500 
501 // chi-squared of the fit
502 inline double TrackBase::chi2() const
503 {
504  return chi2_;
505 }
506 
507 // number of degrees of freedom of the fit
508 inline double TrackBase::ndof() const
509 {
510  return ndof_;
511 }
512 
513 // chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
514 inline double TrackBase::normalizedChi2() const
515 {
516  return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6;
517 }
518 
519 // track electric charge
520 inline int TrackBase::charge() const
521 {
522  return charge_;
523 }
524 
525 // q / p
526 inline double TrackBase::qoverp() const
527 {
528  return charge() / p();
529 }
530 
531 // polar angle
532 inline double TrackBase::theta() const
533 {
534  return momentum_.theta();
535 }
536 
537 // Lambda angle
538 inline double TrackBase::lambda() const
539 {
540  return M_PI_2 - momentum_.theta();
541 }
542 
543 // 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.
544 inline double TrackBase::dxy() const
545 {
546  return (-vx() * py() + vy() * px()) / pt();
547 }
548 
549 // dxy parameter in perigee convention (d0 = -dxy)
550 inline double TrackBase::d0() const
551 {
552  return -dxy();
553 }
554 
555 // 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)
556 inline double TrackBase::dsz() const
557 {
558  return vz() * pt() / p() - (vx() * px() + vy() * py()) / pt() * pz() / p();
559 }
560 
561 // 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.
562 inline double TrackBase::dz() const
563 {
564  return vz() - (vx() * px() + vy() * py()) / pt() * (pz() / pt());
565 }
566 
567 // momentum vector magnitude
568 inline double TrackBase::p() const
569 {
570  return momentum_.R();
571 }
572 
573 // track transverse momentum
574 inline double TrackBase::pt() const
575 {
576  return sqrt(momentum_.Perp2());
577 }
578 
579 // x coordinate of momentum vector
580 inline double TrackBase::px() const
581 {
582  return momentum_.x();
583 }
584 
585 // y coordinate of momentum vector
586 inline double TrackBase::py() const
587 {
588  return momentum_.y();
589 }
590 
591 // z coordinate of momentum vector
592 inline double TrackBase::pz() const
593 {
594  return momentum_.z();
595 }
596 
597 // azimuthal angle of momentum vector
598 inline double TrackBase::phi() const
599 {
600  return momentum_.Phi();
601 }
602 
603 // pseudorapidity of momentum vector
604 inline double TrackBase::eta() const
605 {
606  return momentum_.Eta();
607 }
608 
609 // x coordinate of the reference point on track
610 inline double TrackBase::vx() const
611 {
612  return vertex_.x();
613 }
614 
615 // y coordinate of the reference point on track
616 inline double TrackBase::vy() const
617 {
618  return vertex_.y();
619 }
620 
621 // z coordinate of the reference point on track
622 inline double TrackBase::vz() const
623 {
624  return vertex_.z();
625 }
626 
627 // track momentum vector
629 {
630  return momentum_;
631 }
632 
633 // Reference point on the track
635 {
636  return vertex_;
637 }
638 
639 // reference point on the track. This method is DEPRECATED, please use referencePoint() instead
640 inline const TrackBase::Point & TrackBase::vertex() const
641 {
642  return vertex_;
643 }
644 
645 // dxy parameter with respect to a user-given beamSpot
646 // (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).
647 // This is a good approximation for Tracker tracks.
648 inline double TrackBase::dxy(const Point &myBeamSpot) const
649 {
650  return (-(vx() - myBeamSpot.x()) * py() + (vy() - myBeamSpot.y()) * px()) / pt();
651 }
652 
653 // dxy parameter with respect to the beamSpot taking into account the beamspot slopes
654 // (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).
655 // This is a good approximation for Tracker tracks.
656 inline double TrackBase::dxy(const BeamSpot &theBeamSpot) const
657 {
658  return dxy(theBeamSpot.position(vz()));
659 }
660 
661 // dsz parameter with respect to a user-given beamSpot
662 // (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).
663 // This is a good approximation for Tracker tracks.
664 inline double TrackBase::dsz(const Point &myBeamSpot) const
665 {
666  return (vz() - myBeamSpot.z()) * pt() / p() - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / p();
667 }
668 
669 // dz parameter with respect to a user-given beamSpot
670 // (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).
671 // This is a good approximation for Tracker tracks.
672 inline double TrackBase::dz(const Point &myBeamSpot) const
673 {
674  return (vz() - myBeamSpot.z()) - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / pt();
675 }
676 
677 // Track parameters with one-to-one correspondence to the covariance matrix
679 {
680  return TrackBase::ParameterVector(qoverp(), lambda(), phi(), dxy(), dsz());
681 }
682 
683 // return track covariance matrix
685 {
687  fill(m);
688  return m;
689 }
690 
691 // i-th parameter ( i = 0, ... 4 )
692 inline double TrackBase::parameter(int i) const
693 {
694  return parameters()[i];
695 }
696 
697 // (i,j)-th element of covariance matrix (i, j = 0, ... 4)
698 inline double TrackBase::covariance(int i, int j) const
699 {
700  return covariance_[covIndex(i, j)];
701 }
702 
703 // error on specified element
704 inline double TrackBase::error(int i) const
705 {
706  return sqrt(covariance_[covIndex(i, i)]);
707 }
708 
709 // error on signed transverse curvature
710 inline double TrackBase::qoverpError() const
711 {
712  return error(i_qoverp);
713 }
714 
715 // error on Pt (set to 1000 TeV if charge==0 for safety)
716 inline double TrackBase::ptError() const
717 {
718  return (charge() != 0) ? sqrt(
719  pt() * pt() * p() * p() / charge() / charge() * covariance(i_qoverp, i_qoverp)
720  + 2 * pt() * p() / charge() * pz() * covariance(i_qoverp, i_lambda)
721  + pz() * pz() * covariance(i_lambda, i_lambda)) : 1.e6;
722 }
723 
724 // error on theta
725 inline double TrackBase::thetaError() const
726 {
727  return error(i_lambda);
728 }
729 
730 // error on lambda
731 inline double TrackBase::lambdaError() const
732 {
733  return error(i_lambda);
734 }
735 
736 // error on eta
737 inline double TrackBase::etaError() const
738 {
739  return error(i_lambda) * p() / pt();
740 }
741 
742 // error on phi
743 inline double TrackBase::phiError() const
744 {
745  return error(i_phi);
746 }
747 
748 // error on dxy
749 inline double TrackBase::dxyError() const
750 {
751  return error(i_dxy);
752 }
753 
754 // error on d0
755 inline double TrackBase::d0Error() const
756 {
757  return error(i_dxy);
758 }
759 
760 // error on dsz
761 inline double TrackBase::dszError() const
762 {
763  return error(i_dsz);
764 }
765 
766 // error on dz
767 inline double TrackBase::dzError() const
768 {
769  return error(i_dsz) * p() / pt();
770 }
771 
772 // number of valid hits found
773 inline unsigned short TrackBase::numberOfValidHits() const
774 {
776 }
777 
778 // number of cases where track crossed a layer without getting a hit.
779 inline unsigned short TrackBase::numberOfLostHits() const
780 {
782 }
783 
784 // fraction of valid hits on the track
785 inline double TrackBase::validFraction() const
786 {
791 
792  if ((valid + lost + lostIn + lostOut) == 0) {
793  return -1;
794  }
795 
796  return valid / (double)(valid + lost + lostIn + lostOut);
797 }
798 
799 //Track algorithm
801 {
802  algorithm_ = a;
803  algoMask_.reset();
805 }
806 
808 {
810  algoMask_.set(a);
811 }
812 
813 
814 
815 inline int TrackBase::qualityMask() const
816 {
817  return quality_;
818 }
819 
820 inline void TrackBase::setQualityMask(int qualMask)
821 {
822  quality_ = qualMask;
823 }
824 
825 inline void TrackBase::setNLoops(signed char value)
826 {
827  nLoops_ = value;
828 }
829 
830 inline bool TrackBase::isLooper() const
831 {
832  return (nLoops_ > 0);
833 }
834 
835 inline signed char TrackBase::nLoops() const
836 {
837  return nLoops_;
838 }
839 
840 } // namespace reco
841 
842 #endif
843 
double qoverp() const
q / p
Definition: TrackBase.h:526
double p() const
momentum vector magnitude
Definition: TrackBase.h:568
float chi2_
chi-squared
Definition: TrackBase.h:379
void setQualityMask(int qualMask)
Definition: TrackBase.h:820
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:634
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:755
void setQuality(const TrackQuality)
Definition: TrackBase.h:476
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:485
unsigned int index
index type
Definition: TrackBase.h:95
bool isLooper() const
Definition: TrackBase.h:830
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:443
uint8_t quality_
track quality
Definition: TrackBase.h:404
double validFraction() const
fraction of valid hits on the track
Definition: TrackBase.h:785
uint8_t originalAlgorithm_
track algorithm
Definition: TrackBase.h:400
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:550
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:514
TrackBase()
default constructor
Definition: TrackBase.cc:60
TrackQuality
track quality
Definition: TrackBase.h:139
double theta() const
polar angle
Definition: TrackBase.h:532
double dxyError() const
error on dxy
Definition: TrackBase.h:749
char charge_
electric charge
Definition: TrackBase.h:394
int numberOfValidHits() const
Definition: HitPattern.h:734
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:382
double etaError() const
error on eta
Definition: TrackBase.h:737
#define M_PI_2
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:598
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:779
ErrorD< N >::type type
Definition: Error.h:39
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:74
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:580
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:628
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:556
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:407
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:809
TrackAlgorithm algo() const
Definition: TrackBase.h:450
void setNLoops(signed char value)
Definition: TrackBase.h:825
void setOriginalAlgorithm(const TrackAlgorithm a)
Definition: TrackBase.h:807
bool appendHits(const C &c)
append hit patterns from vector of hit references
Definition: TrackBase.h:438
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:640
double dszError() const
error on dsz
Definition: TrackBase.h:761
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:604
static const std::string qualityNames[]
Definition: TrackBase.h:152
fixed size vector
Definition: Vector.h:31
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:502
std::bitset< algoSize > algoMask_
algo mask, bit set for the algo where it was reconstructed + each algo a track was found overlapping ...
Definition: TrackBase.h:388
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:508
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:684
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:574
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:716
double phiError() const
error on phi
Definition: TrackBase.h:743
int qualityMask() const
Definition: TrackBase.h:815
int j
Definition: DBlmapReader.cc:9
bool appendHitPattern(const TrackingRecHit &hit)
append a single hit to the HitPattern
Definition: TrackBase.h:421
double lambda() const
Lambda angle.
Definition: TrackBase.h:538
const std::complex< double > I
Definition: I.h:8
double error(int i) const
error on specified element
Definition: TrackBase.h:704
bool appendHits(const I &begin, const I &end)
Definition: HitPattern.h:446
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:773
void resetHitPattern()
Sets HitPattern as empty.
Definition: TrackBase.h:426
#define end
Definition: vmac.h:37
math::XYZPoint Point
point in the space
Definition: TrackBase.h:83
void setAlgoMask(AlgoMask a)
Definition: TrackBase.h:333
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:592
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:692
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:710
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:562
double dzError() const
error on dz
Definition: TrackBase.h:767
bool isAlgoInMask(TrackAlgorithm a) const
Definition: TrackBase.h:339
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:622
Definition: DetId.h:18
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:454
std::bitset< algoSize > AlgoMask
algo mask
Definition: TrackBase.h:133
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:114
AlgoMask algoMask() const
Definition: TrackBase.h:335
virtual ~TrackBase()
virtual destructor
Definition: TrackBase.cc:104
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:411
signed char nLoops() const
Definition: TrackBase.h:835
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:461
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:391
static const std::string algoNames[]
Definition: TrackBase.h:136
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:463
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:804
int numberOfValidTrackerHits() const
Definition: HitPattern.h:739
void setAlgorithm(const TrackAlgorithm a)
Track algorithm.
Definition: TrackBase.h:800
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:678
double lambdaError() const
error on lambda
Definition: TrackBase.h:731
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:616
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:376
#define begin
Definition: vmac.h:30
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:373
double a
Definition: hdecay.h:121
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:109
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:126
bool appendHit(const TrackingRecHit &hit)
Definition: HitPattern.cc:168
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:385
int charge() const
track electric charge
Definition: TrackBase.h:520
const Point & position() const
position
Definition: BeamSpot.h:62
unsigned long long algoMaskUL() const
Definition: TrackBase.h:337
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:397
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:544
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:80
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:586
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:610
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
double thetaError() const
error on theta
Definition: TrackBase.h:725