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  // steps used by HI muon regional iterative tracking
139  algoSize = 46
140  };
141 
143  typedef std::bitset<algoSize> AlgoMask;
144 
145 
146  static const std::string algoNames[];
147 
151  loose = 0,
152  tight = 1,
154  confirmed = 3, // means found by more than one iteration
155  goodIterative = 4, // meaningless
158  discarded = 7, // because a better track found. kept in the collection for reference....
160  };
161 
162  static const std::string qualityNames[];
163 
165  TrackBase();
166 
168  TrackBase(double chi2, double ndof, const Point &vertex,
169  const Vector &momentum, int charge, const CovarianceMatrix &cov,
171  signed char nloops = 0);
172 
174  virtual ~TrackBase();
175 
177  double chi2() const;
178 
180  double ndof() const;
181 
183  double normalizedChi2() const;
184 
186  int charge() const;
187 
189  double qoverp() const;
190 
192  double theta() const;
193 
195  double lambda() const;
196 
198  double dxy() const;
199 
201  double d0() const;
202 
204  double dsz() const;
205 
207  double dz() const;
208 
210  double p() const;
211 
213  double pt() const;
214 
216  double px() const;
217 
219  double py() const;
220 
222  double pz() const;
223 
225  double phi() const;
226 
228  double eta() const;
229 
231  double vx() const;
232 
234  double vy() const;
235 
237  double vz() const;
238 
240  const Vector &momentum() const;
241 
243  const Point &referencePoint() const;
244 
246  const Point &vertex() const ;
247  //__attribute__((deprecated("This method is DEPRECATED, please use referencePoint() instead.")));
248 
250  double dxy(const Point &myBeamSpot) const;
251 
253  double dxy(const BeamSpot &theBeamSpot) const;
254 
256  double dsz(const Point &myBeamSpot) const;
257 
259  double dz(const Point &myBeamSpot) const;
260 
262  ParameterVector parameters() const;
263 
265  CovarianceMatrix covariance() const;
266 
268  double parameter(int i) const;
269 
271  double covariance(int i, int j) const;
272 
274  double error(int i) const;
275 
277  double qoverpError() const;
278 
280  double ptError() const;
281 
283  double thetaError() const;
284 
286  double lambdaError() const;
287 
289  double etaError() const;
290 
292  double phiError() const;
293 
295  double dxyError() const;
296 
298  double d0Error() const;
299 
301  double dszError() const;
302 
304  double dzError() const;
305 
307  CovarianceMatrix &fill(CovarianceMatrix &v) const;
308 
310  static index covIndex(index i, index j);
311 
313  const HitPattern &hitPattern() const;
314 
316  unsigned short numberOfValidHits() const;
317 
319  unsigned short numberOfLostHits() const;
320 
322  double validFraction() const;
323 
325  template<typename C>
326  bool appendHits(const C &c);
327 
328  template<typename I>
329  bool appendHits(const I &begin, const I &end);
330 
332  bool appendHitPattern(const TrackingRecHit &hit);
333  bool appendHitPattern(const DetId &id, TrackingRecHit::Type hitType);
334 
336  void resetHitPattern();
337 
339  void setAlgorithm(const TrackAlgorithm a);
340 
342 
343  void setAlgoMask(AlgoMask a) { algoMask_ = a;}
344 
345  AlgoMask algoMask() const { return algoMask_;}
346 #if ( !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__) ) || defined(__ROOTCLING__)
347  unsigned long long algoMaskUL() const { return algoMask().to_ullong();}
348 #endif
349  bool isAlgoInMask(TrackAlgorithm a) const {return algoMask()[a];}
350 
351 
352  TrackAlgorithm algo() const ;
353  TrackAlgorithm originalAlgo() const ;
354 
355 
356  std::string algoName() const;
357 
359 
360  static TrackAlgorithm algoByName(const std::string &name);
361 
363  bool quality(const TrackQuality) const;
364 
365  void setQuality(const TrackQuality);
366 
368 
370 
371  int qualityMask() const;
372 
373  void setQualityMask(int qualMask);
374 
375  void setNLoops(signed char value);
376 
377  bool isLooper() const;
378 
379  signed char nLoops() const;
380 
381 private:
384 
387 
389  float chi2_;
390 
392  Point vertex_;
393 
395  Vector momentum_;
396 
398  std::bitset<algoSize> algoMask_;
399 
401  float ndof_;
402 
404  char charge_;
405 
407  uint8_t algorithm_;
408 
411 
412 
414  uint8_t quality_;
415 
417  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.
418 };
419 
420 // Access the hit pattern, indicating in which Tracker layers the track has hits.
421 inline const HitPattern & TrackBase::hitPattern() const
422 {
423  return hitPattern_;
424 }
425 
427 {
428  return hitPattern_.appendHit(id, hitType);
429 }
430 
432 {
433  return hitPattern_.appendHit(hit);
434 }
435 
437 {
438  hitPattern_.clear();
439 }
440 
441 template<typename I>
442 bool TrackBase::appendHits(const I &begin, const I &end)
443 {
444  return hitPattern_.appendHits(begin, end);
445 }
446 
447 template<typename C>
449 {
450  return setHitPattern(c.begin(), c.end());
451 }
452 
454 {
455  int a = (i <= j ? i : j);
456  int b = (i <= j ? j : i);
457  return b * (b + 1) / 2 + a;
458 }
459 
461 {
462  return (TrackAlgorithm) (algorithm_);
463 }
465 {
467 }
468 
469 
470 
472 
474 {
475  switch (q) {
476  case undefQuality:
477  return quality_ == 0;
478  case goodIterative:
480  default:
481  return (quality_ & (1 << q)) >> q;
482  }
483  return false;
484 }
485 
487 {
488  if (q == undefQuality) {
489  quality_ = 0;
490  } else {
491  quality_ |= (1 << q);
492  }
493 }
494 
496 {
497  if (int(q) < int(qualitySize) && int(q) >= 0) {
498  return qualityNames[int(q)];
499  }
500  return "undefQuality";
501 }
502 
504 {
505  if (int(a) < int(algoSize) && int(a) > 0) {
506  return algoNames[int(a)];
507  }
508  return "undefAlgorithm";
509 }
510 
511 // chi-squared of the fit
512 inline double TrackBase::chi2() const
513 {
514  return chi2_;
515 }
516 
517 // number of degrees of freedom of the fit
518 inline double TrackBase::ndof() const
519 {
520  return ndof_;
521 }
522 
523 // chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
524 inline double TrackBase::normalizedChi2() const
525 {
526  return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6;
527 }
528 
529 // track electric charge
530 inline int TrackBase::charge() const
531 {
532  return charge_;
533 }
534 
535 // q / p
536 inline double TrackBase::qoverp() const
537 {
538  return charge() / p();
539 }
540 
541 // polar angle
542 inline double TrackBase::theta() const
543 {
544  return momentum_.theta();
545 }
546 
547 // Lambda angle
548 inline double TrackBase::lambda() const
549 {
550  return M_PI_2 - momentum_.theta();
551 }
552 
553 // 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.
554 inline double TrackBase::dxy() const
555 {
556  return (-vx() * py() + vy() * px()) / pt();
557 }
558 
559 // dxy parameter in perigee convention (d0 = -dxy)
560 inline double TrackBase::d0() const
561 {
562  return -dxy();
563 }
564 
565 // 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)
566 inline double TrackBase::dsz() const
567 {
568  return vz() * pt() / p() - (vx() * px() + vy() * py()) / pt() * pz() / p();
569 }
570 
571 // 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.
572 inline double TrackBase::dz() const
573 {
574  return vz() - (vx() * px() + vy() * py()) / pt() * (pz() / pt());
575 }
576 
577 // momentum vector magnitude
578 inline double TrackBase::p() const
579 {
580  return momentum_.R();
581 }
582 
583 // track transverse momentum
584 inline double TrackBase::pt() const
585 {
586  return sqrt(momentum_.Perp2());
587 }
588 
589 // x coordinate of momentum vector
590 inline double TrackBase::px() const
591 {
592  return momentum_.x();
593 }
594 
595 // y coordinate of momentum vector
596 inline double TrackBase::py() const
597 {
598  return momentum_.y();
599 }
600 
601 // z coordinate of momentum vector
602 inline double TrackBase::pz() const
603 {
604  return momentum_.z();
605 }
606 
607 // azimuthal angle of momentum vector
608 inline double TrackBase::phi() const
609 {
610  return momentum_.Phi();
611 }
612 
613 // pseudorapidity of momentum vector
614 inline double TrackBase::eta() const
615 {
616  return momentum_.Eta();
617 }
618 
619 // x coordinate of the reference point on track
620 inline double TrackBase::vx() const
621 {
622  return vertex_.x();
623 }
624 
625 // y coordinate of the reference point on track
626 inline double TrackBase::vy() const
627 {
628  return vertex_.y();
629 }
630 
631 // z coordinate of the reference point on track
632 inline double TrackBase::vz() const
633 {
634  return vertex_.z();
635 }
636 
637 // track momentum vector
639 {
640  return momentum_;
641 }
642 
643 // Reference point on the track
645 {
646  return vertex_;
647 }
648 
649 // reference point on the track. This method is DEPRECATED, please use referencePoint() instead
650 inline const TrackBase::Point & TrackBase::vertex() const
651 {
652  return vertex_;
653 }
654 
655 // dxy parameter with respect to a user-given beamSpot
656 // (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).
657 // This is a good approximation for Tracker tracks.
658 inline double TrackBase::dxy(const Point &myBeamSpot) const
659 {
660  return (-(vx() - myBeamSpot.x()) * py() + (vy() - myBeamSpot.y()) * px()) / pt();
661 }
662 
663 // dxy parameter with respect to the beamSpot taking into account the beamspot slopes
664 // (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).
665 // This is a good approximation for Tracker tracks.
666 inline double TrackBase::dxy(const BeamSpot &theBeamSpot) const
667 {
668  return dxy(theBeamSpot.position(vz()));
669 }
670 
671 // dsz parameter with respect to a user-given beamSpot
672 // (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).
673 // This is a good approximation for Tracker tracks.
674 inline double TrackBase::dsz(const Point &myBeamSpot) const
675 {
676  return (vz() - myBeamSpot.z()) * pt() / p() - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / p();
677 }
678 
679 // dz parameter with respect to a user-given beamSpot
680 // (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).
681 // This is a good approximation for Tracker tracks.
682 inline double TrackBase::dz(const Point &myBeamSpot) const
683 {
684  return (vz() - myBeamSpot.z()) - ((vx() - myBeamSpot.x()) * px() + (vy() - myBeamSpot.y()) * py()) / pt() * pz() / pt();
685 }
686 
687 // Track parameters with one-to-one correspondence to the covariance matrix
689 {
690  return TrackBase::ParameterVector(qoverp(), lambda(), phi(), dxy(), dsz());
691 }
692 
693 // return track covariance matrix
695 {
697  fill(m);
698  return m;
699 }
700 
701 // i-th parameter ( i = 0, ... 4 )
702 inline double TrackBase::parameter(int i) const
703 {
704  return parameters()[i];
705 }
706 
707 // (i,j)-th element of covariance matrix (i, j = 0, ... 4)
708 inline double TrackBase::covariance(int i, int j) const
709 {
710  return covariance_[covIndex(i, j)];
711 }
712 
713 // error on specified element
714 inline double TrackBase::error(int i) const
715 {
716  return sqrt(covariance_[covIndex(i, i)]);
717 }
718 
719 // error on signed transverse curvature
720 inline double TrackBase::qoverpError() const
721 {
722  return error(i_qoverp);
723 }
724 
725 // error on Pt (set to 1000 TeV if charge==0 for safety)
726 inline double TrackBase::ptError() const
727 {
728  return (charge() != 0) ? sqrt(
729  pt() * pt() * p() * p() / charge() / charge() * covariance(i_qoverp, i_qoverp)
730  + 2 * pt() * p() / charge() * pz() * covariance(i_qoverp, i_lambda)
731  + pz() * pz() * covariance(i_lambda, i_lambda)) : 1.e6;
732 }
733 
734 // error on theta
735 inline double TrackBase::thetaError() const
736 {
737  return error(i_lambda);
738 }
739 
740 // error on lambda
741 inline double TrackBase::lambdaError() const
742 {
743  return error(i_lambda);
744 }
745 
746 // error on eta
747 inline double TrackBase::etaError() const
748 {
749  return error(i_lambda) * p() / pt();
750 }
751 
752 // error on phi
753 inline double TrackBase::phiError() const
754 {
755  return error(i_phi);
756 }
757 
758 // error on dxy
759 inline double TrackBase::dxyError() const
760 {
761  return error(i_dxy);
762 }
763 
764 // error on d0
765 inline double TrackBase::d0Error() const
766 {
767  return error(i_dxy);
768 }
769 
770 // error on dsz
771 inline double TrackBase::dszError() const
772 {
773  return error(i_dsz);
774 }
775 
776 // error on dz
777 inline double TrackBase::dzError() const
778 {
779  return error(i_dsz) * p() / pt();
780 }
781 
782 // number of valid hits found
783 inline unsigned short TrackBase::numberOfValidHits() const
784 {
786 }
787 
788 // number of cases where track crossed a layer without getting a hit.
789 inline unsigned short TrackBase::numberOfLostHits() const
790 {
792 }
793 
794 // fraction of valid hits on the track
795 inline double TrackBase::validFraction() const
796 {
801 
802  if ((valid + lost + lostIn + lostOut) == 0) {
803  return -1;
804  }
805 
806  return valid / (double)(valid + lost + lostIn + lostOut);
807 }
808 
809 //Track algorithm
811 {
812  algorithm_ = a;
813  algoMask_.reset();
815 }
816 
818 {
820  algoMask_.set(a);
821 }
822 
823 
824 
825 inline int TrackBase::qualityMask() const
826 {
827  return quality_;
828 }
829 
830 inline void TrackBase::setQualityMask(int qualMask)
831 {
832  quality_ = qualMask;
833 }
834 
835 inline void TrackBase::setNLoops(signed char value)
836 {
837  nLoops_ = value;
838 }
839 
840 inline bool TrackBase::isLooper() const
841 {
842  return (nLoops_ > 0);
843 }
844 
845 inline signed char TrackBase::nLoops() const
846 {
847  return nLoops_;
848 }
849 
850 } // namespace reco
851 
852 #endif
853 
double qoverp() const
q / p
Definition: TrackBase.h:536
double p() const
momentum vector magnitude
Definition: TrackBase.h:578
float chi2_
chi-squared
Definition: TrackBase.h:389
void setQualityMask(int qualMask)
Definition: TrackBase.h:830
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:644
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:765
void setQuality(const TrackQuality)
Definition: TrackBase.h:486
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:495
unsigned int index
index type
Definition: TrackBase.h:95
bool isLooper() const
Definition: TrackBase.h:840
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:453
uint8_t quality_
track quality
Definition: TrackBase.h:414
double validFraction() const
fraction of valid hits on the track
Definition: TrackBase.h:795
uint8_t originalAlgorithm_
track algorithm
Definition: TrackBase.h:410
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:560
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:524
TrackBase()
default constructor
Definition: TrackBase.cc:69
TrackQuality
track quality
Definition: TrackBase.h:149
double theta() const
polar angle
Definition: TrackBase.h:542
double dxyError() const
error on dxy
Definition: TrackBase.h:759
char charge_
electric charge
Definition: TrackBase.h:404
int numberOfValidHits() const
Definition: HitPattern.h:734
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:392
double etaError() const
error on eta
Definition: TrackBase.h:747
#define M_PI_2
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:608
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:789
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:590
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:638
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:566
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:417
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:809
TrackAlgorithm algo() const
Definition: TrackBase.h:460
void setNLoops(signed char value)
Definition: TrackBase.h:835
void setOriginalAlgorithm(const TrackAlgorithm a)
Definition: TrackBase.h:817
bool appendHits(const C &c)
append hit patterns from vector of hit references
Definition: TrackBase.h:448
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:650
double dszError() const
error on dsz
Definition: TrackBase.h:771
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:614
static const std::string qualityNames[]
Definition: TrackBase.h:162
fixed size vector
Definition: Vector.h:31
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:512
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:398
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:518
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:694
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:584
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:726
double phiError() const
error on phi
Definition: TrackBase.h:753
int qualityMask() const
Definition: TrackBase.h:825
int j
Definition: DBlmapReader.cc:9
bool appendHitPattern(const TrackingRecHit &hit)
append a single hit to the HitPattern
Definition: TrackBase.h:431
double lambda() const
Lambda angle.
Definition: TrackBase.h:548
const std::complex< double > I
Definition: I.h:8
double error(int i) const
error on specified element
Definition: TrackBase.h:714
bool appendHits(const I &begin, const I &end)
Definition: HitPattern.h:446
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:783
void resetHitPattern()
Sets HitPattern as empty.
Definition: TrackBase.h:436
#define end
Definition: vmac.h:37
math::XYZPoint Point
point in the space
Definition: TrackBase.h:83
void setAlgoMask(AlgoMask a)
Definition: TrackBase.h:343
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:602
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:702
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:720
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:572
double dzError() const
error on dz
Definition: TrackBase.h:777
bool isAlgoInMask(TrackAlgorithm a) const
Definition: TrackBase.h:349
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:632
Definition: DetId.h:18
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:464
std::bitset< algoSize > AlgoMask
algo mask
Definition: TrackBase.h:143
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
AlgoMask algoMask() const
Definition: TrackBase.h:345
virtual ~TrackBase()
virtual destructor
Definition: TrackBase.cc:113
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:421
signed char nLoops() const
Definition: TrackBase.h:845
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:471
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:401
static const std::string algoNames[]
Definition: TrackBase.h:146
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:473
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:810
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:688
double lambdaError() const
error on lambda
Definition: TrackBase.h:741
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:626
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:386
#define begin
Definition: vmac.h:30
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:383
double a
Definition: hdecay.h:121
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:118
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:135
bool appendHit(const TrackingRecHit &hit)
Definition: HitPattern.cc:168
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:395
int charge() const
track electric charge
Definition: TrackBase.h:530
const Point & position() const
position
Definition: BeamSpot.h:62
unsigned long long algoMaskUL() const
Definition: TrackBase.h:347
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:407
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:554
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:80
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:596
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:620
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
double thetaError() const
error on theta
Definition: TrackBase.h:735