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  class TrackBase {
62  public:
64  enum { dimension = 5 };
66  enum { covarianceSize = dimension * ( dimension + 1 ) / 2 };
76  enum { i_qoverp = 0 , i_lambda, i_phi, i_dxy, i_dsz };
78  typedef unsigned int index;
81  iter1=5, iter2=6, iter3=7, iter4=8, iter5=9, iter6=10, iter7=11, iter8=12, iter9=13,iter10=14,
87  beamhalo=28,
88  gsf=29,
89  algoSize=30 };
90  static const std::string algoNames[];
91 
94  static const std::string qualityNames[];
95 
97  TrackBase();
99  TrackBase( double chi2, double ndof, const Point & referencePoint,
100  const Vector & momentum, int charge, const CovarianceMatrix &,
103  virtual ~TrackBase();
105  double chi2() const { return chi2_; }
107  double ndof() const { return ndof_; }
109  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
111  int charge() const { return charge_; }
113  double qoverp() const { return charge() / p(); }
115  double theta() const { return momentum_.theta(); }
117  double lambda() const { return M_PI/2 - momentum_.theta(); }
119  double dxy() const { return ( - vx() * py() + vy() * px() ) / pt(); }
121  double d0() const { return - dxy(); }
123  double dsz() const { return vz()*pt()/p() - (vx()*px()+vy()*py())/pt() * pz()/p(); }
125  double dz() const { return vz() - (vx()*px()+vy()*py())/pt() * (pz()/pt()); }
127  double p() const { return momentum_.R(); }
129  double pt() const { return sqrt( momentum_.Perp2() ); }
131  double px() const { return momentum_.x(); }
133  double py() const { return momentum_.y(); }
135  double pz() const { return momentum_.z(); }
137  double phi() const { return momentum_.Phi(); }
139  double eta() const { return momentum_.Eta(); }
141  double vx() const { return vertex_.x(); }
143  double vy() const { return vertex_.y(); }
145  double vz() const { return vertex_.z(); }
146 
148  const Vector & momentum() const { return momentum_; }
149 
151  const Point & referencePoint() const { return vertex_; }
152 
154  const Point & vertex() const { return vertex_; }
155 
157  double dxy(const Point& myBeamSpot) const {
158  return ( - (vx()-myBeamSpot.x()) * py() + (vy()-myBeamSpot.y()) * px() ) / pt();
159  }
160 
162  double dxy(const BeamSpot& theBeamSpot) const {
163  return dxy(theBeamSpot.position(vz()));
164  }
165 
167  double dsz(const Point& myBeamSpot) const {
168  return (vz()-myBeamSpot.z())*pt()/p() - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() *pz()/p();
169  }
171  double dz(const Point& myBeamSpot) const {
172  return (vz()-myBeamSpot.z()) - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() * pz()/pt();
173  }
174 
176  ParameterVector parameters() const {
177  return ParameterVector(qoverp(),lambda(),phi(),dxy(),dsz());
178  }
180  CovarianceMatrix covariance() const { CovarianceMatrix m; fill( m ); return m; }
181 
183  double parameter(int i) const { return parameters()[i]; }
185  double covariance( int i, int j ) const { return covariance_[ covIndex( i, j ) ]; }
187  double error( int i ) const { return sqrt( covariance_[ covIndex( i, i ) ] ); }
188 
190  double qoverpError() const { return error( i_qoverp ); }
192  double ptError() const {
193  return (charge()!=0) ? sqrt(
194  pt()*pt()*p()*p()/charge()/charge() * covariance(i_qoverp,i_qoverp)
195  + 2*pt()*p()/charge()*pz() * covariance(i_qoverp,i_lambda)
196  + pz()*pz() * covariance(i_lambda,i_lambda) ) : 1.e6;
197  }
199  double thetaError() const { return error( i_lambda ); }
201  double lambdaError() const { return error( i_lambda ); }
203  double etaError() const { return error( i_lambda ) * p()/pt(); }
205  double phiError() const { return error( i_phi ); }
207  double dxyError() const { return error( i_dxy ); }
209  double d0Error() const { return error( i_dxy ); }
211  double dszError() const { return error( i_dsz ); }
213  double dzError() const { return error( i_dsz ) * p()/pt(); }
214 
216  CovarianceMatrix & fill( CovarianceMatrix & v ) const;
218  static index covIndex( index i, index j );
219 
221  const HitPattern & hitPattern() const { return hitPattern_; }
226 
227 
228 
230  unsigned short numberOfValidHits() const { return hitPattern_.numberOfValidHits(); }
232  unsigned short numberOfLostHits() const { return hitPattern_.numberOfLostHits(); }
234  double validFraction() const {
239  if ((valid+lost+lostIn+lostOut)==0) return -1;
240  return valid/(1.0*(valid+lost+lostIn+lostOut));
241  }
243  template<typename C>
244  void setHitPattern( const C & c ) { hitPattern_.set( c.begin(), c.end() ); }
245  template<typename C>
246  void setTrackerExpectedHitsInner( const C & c ) { trackerExpectedHitsInner_.set( c.begin(), c.end() ); }
247  template<typename C>
248  void setTrackerExpectedHitsOuter( const C & c ) { trackerExpectedHitsOuter_.set( c.begin(), c.end() ); }
249 
250  template<typename I>
251  void setHitPattern( const I & begin, const I & end ) { hitPattern_.set( begin, end ); }
252  template<typename I>
253  void setTrackerExpectedHitsInner( const I & begin, const I & end ) { trackerExpectedHitsInner_.set( begin, end ); }
254  template<typename I>
255  void setTrackerExpectedHitsOuter( const I & begin, const I & end ) { trackerExpectedHitsOuter_.set( begin, end ); }
256 
258  void setHitPattern( const TrackingRecHit & hit, size_t i ) { hitPattern_.set( hit, i ); }
262 
264  void setHitPattern( const HitPattern& hitP) {hitPattern_ = hitP;}
267 
269 
271  void setAlgorithm(const TrackAlgorithm a, bool set=true) { if (set) algorithm_=a; else algorithm_=undefAlgorithm;}
272  TrackAlgorithm algo() const ;
273  std::string algoName() const;
275  static TrackAlgorithm algoByName(const std::string &name);
276 
278  bool quality(const TrackQuality ) const;
279  void setQuality(const TrackQuality, bool set=true);
282 
283 
284  int qualityMask() const { return quality_; }
285  void setQualityMask(int qualMask) {quality_ = qualMask;}
286 
287 
288  void setNLoops(signed char value) { nLoops_=value;}
289 
290  bool isLooper() const { return (nLoops_>0);}
291  signed char nLoops() const {return nLoops_;}
292 
293  private:
300 
304  float chi2_;
306  Point vertex_;
308  Vector momentum_;
309 
311  float ndof_;
312 
314  char charge_;
316  uint8_t algorithm_;
318  uint8_t quality_;
319 
321  signed char nLoops_;
322 
323  };
324 
326  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
327  return b * ( b + 1 ) / 2 + a;
328  }
329 
331  return (TrackAlgorithm) algorithm_;
332  }
333 
335  // I'd like to do:
336  // return TrackBase::algoName(algorithm_);
337  // but I cannot define a const static function. Why???
338 
339  switch(algorithm_)
340  {
341  case undefAlgorithm: return "undefAlgorithm";
342  case ctf: return "ctf";
343  case rs: return "rs";
344  case cosmics: return "cosmics";
345  case beamhalo: return "beamhalo";
346  case iter0: return "iter0";
347  case iter1: return "iter1";
348  case iter2: return "iter2";
349  case iter3: return "iter3";
350  case iter4: return "iter4";
351  case iter5: return "iter5";
352  case iter6: return "iter6";
353  case iter7: return "iter7";
354  case iter8: return "iter8";
355  case iter9: return "iter9";
356  case iter10: return "iter10";
357  case outInEcalSeededConv: return "outInEcalSeededConv";
358  case inOutEcalSeededConv: return "inOutEcalSeededConv";
359  case nuclInter: return "nuclInter";
360  case standAloneMuon: return "standAloneMuon";
361  case globalMuon: return "globalMuon";
362  case cosmicStandAloneMuon: return "cosmicStandAloneMuon";
363  case cosmicGlobalMuon: return "cosmicGlobalMuon";
364  case iter1LargeD0: return "iter1LargeD0";
365  case iter2LargeD0: return "iter2LargeD0";
366  case iter3LargeD0: return "iter3LargeD0";
367  case iter4LargeD0: return "iter4LargeD0";
368  case iter5LargeD0: return "iter5LargeD0";
369  case bTagGhostTracks: return "bTagGhostTracks";
370  case gsf: return "gsf";
371  }
372  return "undefAlgorithm";
373  }
374 
375  inline bool TrackBase::quality( const TrackBase::TrackQuality q) const{
376  switch(q){
377  case undefQuality: return (quality_==0);
378  case goodIterative:
379  {
380  return ( ((quality_ & (1<<TrackBase::confirmed))>>TrackBase::confirmed) ||
382  }
383  default:
384  {
385  return (quality_ & (1<<q))>>q;
386  }
387  }
388  return false;
389  }
390 
391  inline void TrackBase::setQuality( const TrackBase::TrackQuality q, bool set){
392  switch(q){
393  case undefQuality: quality_=0;
394  default:
395  {
396  if (set)//regular OR if setting value to true
397  quality_= quality_ | (1<<q) ;
398  else // doing "half-XOR" if unsetting value
399  quality_= quality_ & (~(1<<q));
400  }
401  }
402  }
403 
405  if(int(q) < int(qualitySize) && int(q)>=0) return qualityNames[int(q)];
406  return "undefQuality";
407  }
408 
410  if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
411  return "undefAlgorithm";
412  }
413 
414 
415 
416 }
417 
418 #endif
double qoverp() const
q/p
Definition: TrackBase.h:113
double p() const
momentum vector magnitude
Definition: TrackBase.h:127
float chi2_
chi-squared
Definition: TrackBase.h:304
void setQualityMask(int qualMask)
Definition: TrackBase.h:285
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:209
double dxy(const Point &myBeamSpot) const
dxy parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted a...
Definition: TrackBase.h:157
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:404
unsigned int index
index type
Definition: TrackBase.h:78
bool isLooper() const
Definition: TrackBase.h:290
void setTrackerExpectedHitsInner(const C &c)
Definition: TrackBase.h:246
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:325
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:148
uint8_t quality_
track quality
Definition: TrackBase.h:318
double validFraction() const
fraction of valid hits on the track
Definition: TrackBase.h:234
double d0() const
dxy parameter in perigee convention (d0 = - dxy)
Definition: TrackBase.h:121
void setTrackerExpectedHitsOuter(const I &begin, const I &end)
Definition: TrackBase.h:255
void appendHit(const TrackingRecHit &hit)
Definition: HitPattern.cc:86
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:109
TrackBase()
default constructor
Definition: TrackBase.cc:20
double dxy(const BeamSpot &theBeamSpot) const
dxy parameter with respect to the beamSpot taking into account the beamspot slopes WARNING: this quan...
Definition: TrackBase.h:162
int numberOfLostTrackerHits() const
Definition: HitPattern.h:650
TrackQuality
track quality
Definition: TrackBase.h:93
double theta() const
polar angle
Definition: TrackBase.h:115
double dxyError() const
error on dxy
Definition: TrackBase.h:207
char charge_
electric charge
Definition: TrackBase.h:314
int numberOfValidHits() const
Definition: HitPattern.h:589
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:306
void setTrackerExpectedHitsInner(const I &begin, const I &end)
Definition: TrackBase.h:253
int numberOfLostHits() const
Definition: HitPattern.h:646
void setHitPattern(const TrackingRecHit &hit, size_t i)
set hit pattern for specified hit
Definition: TrackBase.h:258
void setHitPattern(const HitPattern &hitP)
set hitPattern from pre-defined hitPattern
Definition: TrackBase.h:264
double etaError() const
error on eta
Definition: TrackBase.h:203
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:137
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:232
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:297
ErrorD< N >::type type
Definition: Error.h:39
void set(const I &begin, const I &end)
Definition: HitPattern.h:152
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:68
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:131
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:151
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:123
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:321
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:180
TrackAlgorithm
track algorithm
Definition: TrackBase.h:80
TrackAlgorithm algo() const
Definition: TrackBase.h:330
void setTrackerExpectedHitsInner(const TrackingRecHit &hit, size_t i)
Definition: TrackBase.h:260
void setNLoops(signed char value)
Definition: TrackBase.h:288
void setHitPattern(const I &begin, const I &end)
Definition: TrackBase.h:251
double dszError() const
error on dsz
Definition: TrackBase.h:211
const HitPattern & trackerExpectedHitsOuter() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers after the last...
Definition: TrackBase.h:225
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
void setQuality(const TrackQuality, bool set=true)
Definition: TrackBase.h:391
static const std::string qualityNames[]
Definition: TrackBase.h:94
fixed size vector
Definition: Vector.h:31
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:105
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:176
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:107
T sqrt(T t)
Definition: SSEVec.h:48
void setTrackerExpectedHitsOuter(const C &c)
Definition: TrackBase.h:248
double pt() const
track transverse momentum
Definition: TrackBase.h:129
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:192
double phiError() const
error on phi
Definition: TrackBase.h:205
void setTrackerExpectedHitsOuter(const TrackingRecHit &hit, size_t i)
Definition: TrackBase.h:261
int qualityMask() const
Definition: TrackBase.h:284
int j
Definition: DBlmapReader.cc:9
double lambda() const
Lambda angle.
Definition: TrackBase.h:117
const std::complex< double > I
Definition: I.h:8
double error(int i) const
error on specified element
Definition: TrackBase.h:187
const HitPattern & trackerExpectedHitsInner() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the fir...
Definition: TrackBase.h:223
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:230
#define end
Definition: vmac.h:37
math::XYZPoint Point
point in the space
Definition: TrackBase.h:74
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
#define M_PI
void setTrackerExpectedHitsOuter(const HitPattern &hitP)
Definition: TrackBase.h:266
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:135
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:183
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:190
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:125
double dzError() const
error on dz
Definition: TrackBase.h:213
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:154
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:145
void setHitPattern(const C &c)
set hit patterns from vector of hit references
Definition: TrackBase.h:244
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
virtual ~TrackBase()
virtual destructor
Definition: TrackBase.cc:39
signed char nLoops() const
Definition: TrackBase.h:291
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:334
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:311
static const std::string algoNames[]
Definition: TrackBase.h:90
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:375
double dsz(const Point &myBeamSpot) const
dsz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted a...
Definition: TrackBase.h:167
int numberOfValidTrackerHits() const
Definition: HitPattern.h:593
double lambdaError() const
error on lambda
Definition: TrackBase.h:201
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:143
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:302
#define begin
Definition: vmac.h:30
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:295
void setAlgorithm(const TrackAlgorithm a, bool set=true)
position index
Definition: TrackBase.h:271
double a
Definition: hdecay.h:121
double covariance(int i, int j) const
(i,j)-th element of covarianve matrix ( i, j = 0, ... 4 )
Definition: TrackBase.h:185
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:42
void appendHitPattern(const TrackingRecHit &hit)
Definition: TrackBase.h:259
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:299
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:55
double dz(const Point &myBeamSpot) const
dz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as...
Definition: TrackBase.h:171
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:308
int charge() const
track electric charge
Definition: TrackBase.h:111
const Point & position() const
position
Definition: BeamSpot.h:62
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:316
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:119
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:72
void setTrackerExpectedHitsInner(const HitPattern &hitP)
Definition: TrackBase.h:265
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:133
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:141
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:70
double thetaError() const
error on theta
Definition: TrackBase.h:199