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 
59 
60 namespace reco {
61 
62  class TrackBase {
63  public:
65  enum { dimension = 5 };
67  enum { covarianceSize = dimension * ( dimension + 1 ) / 2 };
77  enum { i_qoverp = 0 , i_lambda, i_phi, i_dxy, i_dsz };
79  typedef unsigned int index;
82  iter1=5, iter2=6, iter3=7, iter4=8, iter5=9, iter6=10, iter7=11, iter8=12, iter9=13,iter10=14,
88  beamhalo=28,
89  gsf=29,
90  algoSize=30 };
91  static const std::string algoNames[];
92 
95  static const std::string qualityNames[];
96 
98  TrackBase();
100  TrackBase( double chi2, double ndof, const Point & referencePoint,
101  const Vector & momentum, int charge, const CovarianceMatrix &,
104  ~TrackBase();
106  double chi2() const { return chi2_; }
108  double ndof() const { return ndof_; }
110  double normalizedChi2() const { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
112  int charge() const { return charge_; }
114  double qoverp() const { return charge() / p(); }
116  double theta() const { return momentum_.theta(); }
118  double lambda() const { return M_PI/2 - momentum_.theta(); }
120  double dxy() const { return ( - vx() * py() + vy() * px() ) / pt(); }
122  double d0() const { return - dxy(); }
124  double dsz() const { return vz()*pt()/p() - (vx()*px()+vy()*py())/pt() * pz()/p(); }
126  double dz() const { return vz() - (vx()*px()+vy()*py())/pt() * (pz()/pt()); }
128  double p() const { return momentum_.R(); }
130  double pt() const { return sqrt( momentum_.Perp2() ); }
132  double px() const { return momentum_.x(); }
134  double py() const { return momentum_.y(); }
136  double pz() const { return momentum_.z(); }
138  double phi() const { return momentum_.Phi(); }
140  double eta() const { return momentum_.Eta(); }
142  double vx() const { return vertex_.x(); }
144  double vy() const { return vertex_.y(); }
146  double vz() const { return vertex_.z(); }
147 
149  const Vector & momentum() const { return momentum_; }
150 
152  const Point & referencePoint() const { return vertex_; }
153 
155  const Point & vertex() const { return vertex_; }
156 
158  double dxy(const Point& myBeamSpot) const {
159  return ( - (vx()-myBeamSpot.x()) * py() + (vy()-myBeamSpot.y()) * px() ) / pt();
160  }
162  double dsz(const Point& myBeamSpot) const {
163  return (vz()-myBeamSpot.z())*pt()/p() - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() *pz()/p();
164  }
166  double dz(const Point& myBeamSpot) const {
167  return (vz()-myBeamSpot.z()) - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() * pz()/pt();
168  }
169 
171  ParameterVector parameters() const {
172  return ParameterVector(qoverp(),lambda(),phi(),dxy(),dsz());
173  }
175  CovarianceMatrix covariance() const { CovarianceMatrix m; fill( m ); return m; }
176 
178  double parameter(int i) const { return parameters()[i]; }
180  double covariance( int i, int j ) const { return covariance_[ covIndex( i, j ) ]; }
182  double error( int i ) const { return sqrt( covariance_[ covIndex( i, i ) ] ); }
183 
185  double qoverpError() const { return error( i_qoverp ); }
187  double ptError() const {
188  return (charge()!=0) ? sqrt(
189  pt()*pt()*p()*p()/charge()/charge() * covariance(i_qoverp,i_qoverp)
190  + 2*pt()*p()/charge()*pz() * covariance(i_qoverp,i_lambda)
191  + pz()*pz() * covariance(i_lambda,i_lambda) ) : 1.e6;
192  }
194  double thetaError() const { return error( i_lambda ); }
196  double lambdaError() const { return error( i_lambda ); }
198  double etaError() const { return error( i_lambda ) * p()/pt(); }
200  double phiError() const { return error( i_phi ); }
202  double dxyError() const { return error( i_dxy ); }
204  double d0Error() const { return error( i_dxy ); }
206  double dszError() const { return error( i_dsz ); }
208  double dzError() const { return error( i_dsz ) * p()/pt(); }
209 
211  CovarianceMatrix & fill( CovarianceMatrix & v ) const;
213  static index covIndex( index i, index j );
214 
216  const HitPattern & hitPattern() const { return hitPattern_; }
221 
222 
223 
225  unsigned short numberOfValidHits() const { return hitPattern_.numberOfValidHits(); }
227  unsigned short numberOfLostHits() const { return hitPattern_.numberOfLostHits(); }
229  template<typename C>
230  void setHitPattern( const C & c ) { hitPattern_.set( c.begin(), c.end() ); }
231  template<typename C>
232  void setTrackerExpectedHitsInner( const C & c ) { trackerExpectedHitsInner_.set( c.begin(), c.end() ); }
233  template<typename C>
234  void setTrackerExpectedHitsOuter( const C & c ) { trackerExpectedHitsOuter_.set( c.begin(), c.end() ); }
235 
236  template<typename I>
237  void setHitPattern( const I & begin, const I & end ) { hitPattern_.set( begin, end ); }
238  template<typename I>
239  void setTrackerExpectedHitsInner( const I & begin, const I & end ) { trackerExpectedHitsInner_.set( begin, end ); }
240  template<typename I>
241  void setTrackerExpectedHitsOuter( const I & begin, const I & end ) { trackerExpectedHitsOuter_.set( begin, end ); }
242 
244  void setHitPattern( const TrackingRecHit & hit, size_t i ) { hitPattern_.set( hit, i ); }
247 
249  void setHitPattern( const HitPattern& hitP) {hitPattern_ = hitP;}
252 
254 
256  void setAlgorithm(const TrackAlgorithm a, bool set=true) { if (set) algorithm_=a; else algorithm_=undefAlgorithm;}
257  TrackAlgorithm algo() const ;
258  std::string algoName() const;
259  static std::string algoName(TrackAlgorithm);
260  static TrackAlgorithm algoByName(const std::string &name);
261 
263  bool quality(const TrackQuality ) const;
264  void setQuality(const TrackQuality, bool set=true);
265  static std::string qualityName(TrackQuality);
266  static TrackQuality qualityByName(const std::string &name);
267 
268 
269  int qualityMask() const { return quality_; }
270  void setQualityMask(int qualMask) {quality_ = qualMask;}
271  private:
273  float chi2_;
275  float ndof_;
277  Point vertex_;
279  Vector momentum_;
281  char charge_;
290 
292  uint8_t algorithm_;
294  uint8_t quality_;
295 
296  };
297 
299  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
300  return b * ( b + 1 ) / 2 + a;
301  }
302 
304  return (TrackAlgorithm) algorithm_;
305  }
306 
307  inline std::string TrackBase::algoName() const{
308  // I'd like to do:
309  // return TrackBase::algoName(algorithm_);
310  // but I cannot define a const static function. Why???
311 
312  switch(algorithm_)
313  {
314  case undefAlgorithm: return "undefAlgorithm";
315  case ctf: return "ctf";
316  case rs: return "rs";
317  case cosmics: return "cosmics";
318  case beamhalo: return "beamhalo";
319  case iter0: return "iter0";
320  case iter1: return "iter1";
321  case iter2: return "iter2";
322  case iter3: return "iter3";
323  case iter4: return "iter4";
324  case iter5: return "iter5";
325  case iter6: return "iter6";
326  case iter7: return "iter7";
327  case iter8: return "iter8";
328  case iter9: return "iter9";
329  case iter10: return "iter10";
330  case outInEcalSeededConv: return "outInEcalSeededConv";
331  case inOutEcalSeededConv: return "inOutEcalSeededConv";
332  case nuclInter: return "nuclInter";
333  case standAloneMuon: return "standAloneMuon";
334  case globalMuon: return "globalMuon";
335  case cosmicStandAloneMuon: return "cosmicStandAloneMuon";
336  case cosmicGlobalMuon: return "cosmicGlobalMuon";
337  case iter1LargeD0: return "iter1LargeD0";
338  case iter2LargeD0: return "iter2LargeD0";
339  case iter3LargeD0: return "iter3LargeD0";
340  case iter4LargeD0: return "iter4LargeD0";
341  case iter5LargeD0: return "iter5LargeD0";
342  case bTagGhostTracks: return "bTagGhostTracks";
343  case gsf: return "gsf";
344  }
345  return "undefAlgorithm";
346  }
347 
348  inline bool TrackBase::quality( const TrackBase::TrackQuality q) const{
349  switch(q){
350  case undefQuality: return (quality_==0);
351  case goodIterative:
352  {
353  return ( ((quality_ & (1<<TrackBase::confirmed))>>TrackBase::confirmed) ||
355  }
356  default:
357  {
358  return (quality_ & (1<<q))>>q;
359  }
360  }
361  return false;
362  }
363 
364  inline void TrackBase::setQuality( const TrackBase::TrackQuality q, bool set){
365  switch(q){
366  case undefQuality: quality_=0;
367  default:
368  {
369  if (set)//regular OR if setting value to true
370  quality_= quality_ | (1<<q) ;
371  else // doing "half-XOR" if unsetting value
372  quality_= quality_ & (~(1<<q));
373  }
374  }
375  }
376 
377  inline std::string TrackBase::qualityName(TrackQuality q) {
378  if(int(q) < int(qualitySize) && int(q)>=0) return qualityNames[int(q)];
379  return "undefQuality";
380  }
381 
382  inline std::string TrackBase::algoName(TrackAlgorithm a){
383  if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
384  return "undefAlgorithm";
385  }
386 
387 
388 
389 }
390 
391 #endif
double qoverp() const
q/p
Definition: TrackBase.h:114
double p() const
momentum vector magnitude
Definition: TrackBase.h:128
float chi2_
chi-squared
Definition: TrackBase.h:273
void setQualityMask(int qualMask)
Definition: TrackBase.h:270
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:204
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:158
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:377
unsigned int index
index type
Definition: TrackBase.h:79
int numberOfValidHits() const
Definition: HitPattern.cc:321
void setTrackerExpectedHitsInner(const C &c)
Definition: TrackBase.h:232
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:298
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:149
uint8_t quality_
track quality
Definition: TrackBase.h:294
double d0() const
dxy parameter in perigee convention (d0 = - dxy)
Definition: TrackBase.h:122
void setTrackerExpectedHitsOuter(const I &begin, const I &end)
Definition: TrackBase.h:241
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:110
TrackBase()
default constructor
Definition: TrackBase.cc:20
TrackQuality
track quality
Definition: TrackBase.h:94
double theta() const
polar angle
Definition: TrackBase.h:116
double dxyError() const
error on dxy
Definition: TrackBase.h:202
char charge_
electric charge
Definition: TrackBase.h:281
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:277
void setTrackerExpectedHitsInner(const I &begin, const I &end)
Definition: TrackBase.h:239
void setHitPattern(const TrackingRecHit &hit, size_t i)
set hit pattern for specified hit
Definition: TrackBase.h:244
void setHitPattern(const HitPattern &hitP)
set hitPattern from pre-defined hitPattern
Definition: TrackBase.h:249
double etaError() const
error on eta
Definition: TrackBase.h:198
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:138
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:227
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:287
ErrorD< N >::type type
Definition: Error.h:30
void set(const I &begin, const I &end)
Definition: HitPattern.h:138
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:69
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:132
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:152
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:124
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:175
TrackAlgorithm
track algorithm
Definition: TrackBase.h:81
TrackAlgorithm algo() const
Definition: TrackBase.h:303
void setTrackerExpectedHitsInner(const TrackingRecHit &hit, size_t i)
Definition: TrackBase.h:245
void setHitPattern(const I &begin, const I &end)
Definition: TrackBase.h:237
double dszError() const
error on dsz
Definition: TrackBase.h:206
const HitPattern & trackerExpectedHitsOuter() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers after the last...
Definition: TrackBase.h:220
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:140
void setQuality(const TrackQuality, bool set=true)
Definition: TrackBase.h:364
static const std::string qualityNames[]
Definition: TrackBase.h:95
fixed size vector
Definition: Vector.h:23
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:106
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:171
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:108
T sqrt(T t)
Definition: SSEVec.h:28
void setTrackerExpectedHitsOuter(const C &c)
Definition: TrackBase.h:234
double pt() const
track transverse momentum
Definition: TrackBase.h:130
int numberOfLostHits() const
Definition: HitPattern.cc:502
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:187
double phiError() const
error on phi
Definition: TrackBase.h:200
void setTrackerExpectedHitsOuter(const TrackingRecHit &hit, size_t i)
Definition: TrackBase.h:246
int qualityMask() const
Definition: TrackBase.h:269
int j
Definition: DBlmapReader.cc:9
double lambda() const
Lambda angle.
Definition: TrackBase.h:118
const std::complex< double > I
Definition: I.h:8
double error(int i) const
error on specified element
Definition: TrackBase.h:182
const HitPattern & trackerExpectedHitsInner() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the fir...
Definition: TrackBase.h:218
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:225
#define end
Definition: vmac.h:38
math::XYZPoint Point
point in the space
Definition: TrackBase.h:75
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:216
void setTrackerExpectedHitsOuter(const HitPattern &hitP)
Definition: TrackBase.h:251
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:136
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:178
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:185
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:126
double dzError() const
error on dz
Definition: TrackBase.h:208
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:155
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:146
void setHitPattern(const C &c)
set hit patterns from vector of hit references
Definition: TrackBase.h:230
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
#define M_PI
Definition: BFit3D.cc:3
~TrackBase()
virtual destructor
Definition: TrackBase.cc:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
std::string algoName() const
Definition: TrackBase.h:307
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:26
double b
Definition: hdecay.h:120
float ndof_
number of degrees of freedom
Definition: TrackBase.h:275
static const std::string algoNames[]
Definition: TrackBase.h:91
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:348
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:162
double lambdaError() const
error on lambda
Definition: TrackBase.h:196
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:144
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:283
#define begin
Definition: vmac.h:31
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:285
void setAlgorithm(const TrackAlgorithm a, bool set=true)
position index
Definition: TrackBase.h:256
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:180
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:42
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:289
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:166
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:279
int charge() const
track electric charge
Definition: TrackBase.h:112
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:292
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:120
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:73
void setTrackerExpectedHitsInner(const HitPattern &hitP)
Definition: TrackBase.h:250
mathSSE::Vec4< T > v
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:134
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:142
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:71
double thetaError() const
error on theta
Definition: TrackBase.h:194