CMS 3D CMS Logo

TrackBase.cc
Go to the documentation of this file.
1 #include "Rtypes.h"
4 #include <algorithm>
5 
6 using namespace reco;
7 
8 // To be kept in synch with the enumerator definitions in TrackBase.h file
10  "undefAlgorithm",
11  "ctf",
12  "duplicateMerge",
13  "cosmics",
14  "initialStep",
15  "lowPtTripletStep",
16  "pixelPairStep",
17  "detachedTripletStep",
18  "mixedTripletStep",
19  "pixelLessStep",
20  "tobTecStep",
21  "jetCoreRegionalStep",
22  "conversionStep",
23  "muonSeededStepInOut",
24  "muonSeededStepOutIn",
25  "outInEcalSeededConv",
26  "inOutEcalSeededConv",
27  "nuclInter",
28  "standAloneMuon",
29  "globalMuon",
30  "cosmicStandAloneMuon",
31  "cosmicGlobalMuon",
32  "highPtTripletStep",
33  "lowPtQuadStep",
34  "detachedQuadStep",
35  "reservedForUpgrades1",
36  "reservedForUpgrades2",
37  "bTagGhostTracks",
38  "beamhalo" ,
39  "gsf",
40  "hltPixel",
41  "hltIter0",
42  "hltIter1",
43  "hltIter2",
44  "hltIter3",
45  "hltIter4",
46  "hltIterX",
47  "hiRegitMuInitialStep",
48  "hiRegitMuLowPtTripletStep",
49  "hiRegitMuPixelPairStep",
50  "hiRegitMuDetachedTripletStep",
51  "hiRegitMuMixedTripletStep",
52  "hiRegitMuPixelLessStep",
53  "hiRegitMuTobTecStep",
54  "hiRegitMuMuonSeededStepInOut",
55  "hiRegitMuMuonSeededStepOutIn"
56 };
57 
59  "loose",
60  "tight",
61  "highPurity",
62  "confirmed",
63  "goodIterative",
64  "looseSetWithPV",
65  "highPuritySetWithPV",
66  "discarded"
67 };
68 
70  covt0t0_(-1.f),
71  covbetabeta_(-1.f),
72  chi2_(0),
73  vertex_(0, 0, 0),
74  t0_(0),
75  momentum_(0, 0, 0),
76  beta_(0),
77  ndof_(0),
78  charge_(0),
79  algorithm_(undefAlgorithm),
80  originalAlgorithm_(undefAlgorithm),
81  quality_(0),
82  nLoops_(0),
83  stopReason_(0)
84 {
85  algoMask_.set(algorithm_);
86  index idx = 0;
87  for (index i = 0; i < dimension; ++i) {
88  for (index j = 0; j <= i; ++j) {
89  covariance_[idx++] = 0;
90  }
91  }
92 }
93 
94 TrackBase::TrackBase(double chi2, double ndof, const Point &vertex, const Vector &momentum,
96  TrackQuality quality, signed char nloops, uint8_t stopReason,
97  float t0, float beta, float covt0t0, float covbetabeta):
98  covt0t0_(covt0t0),
99  covbetabeta_(covbetabeta),
100  chi2_(chi2),
101  vertex_(vertex),
102  t0_(t0),
103  momentum_(momentum),
104  beta_(beta),
105  ndof_(ndof),
106  charge_(charge),
107  algorithm_(algorithm),
108  originalAlgorithm_(algorithm),
109  quality_(0),
110  nLoops_(nloops),
111  stopReason_(stopReason)
112 {
113  algoMask_.set(algorithm_);
114 
115  index idx = 0;
116  for (index i = 0; i < dimension; ++i) {
117  for (index j = 0; j <= i; ++j) {
118  covariance_[idx++] = cov(i, j);
119  }
120  }
121  setQuality(quality);
122 }
123 
125 {
126  ;
127 }
128 
130 {
131  return fillCovariance(v, covariance_);
132 }
133 
135 {
137  int index = std::find(qualityNames, qualityNames + size, name) - qualityNames;
138  if (index == size) {
139  return undefQuality; // better this or throw() ?
140  }
141 
142  // cast
143  return TrackQuality(index);
144 }
145 
147 {
149  int index = std::find(algoNames, algoNames + size, name) - algoNames;
150  if (index == size) {
151  return undefAlgorithm; // better this or throw() ?
152  }
153 
154  // cast
155  return TrackAlgorithm(index);
156 }
157 
158 double TrackBase::dxyError(Point const &vtx, math::Error<3>::type const &vertexCov) const {
159  // Gradient of TrackBase::dxy(const Point &myBeamSpot) with respect to track parameters. Using unrolled expressions to avoid calling for higher dimension matrices
160  // ( 0, 0, x_vert * cos(phi) + y_vert * sin(phi), 1, 0 )
161  // Gradient with respect to point parameters
162  // ( sin(phi), -cos(phi))
163  // Propagate covariance assuming cross-terms of the covariance between track and vertex parameters are 0
164  return std::sqrt((vtx.x() * px() + vtx.y() * py()) * (vtx.x() * px() + vtx.y() * py()) / (pt() * pt()) *
166  2 * (vtx.x() * px() + vtx.y() * py()) / pt() * covariance(i_phi, i_dxy) + covariance(i_dxy, i_dxy) +
167  py() * py() / (pt() * pt()) * vertexCov(0, 0) - 2 * py() * px() / (pt() * pt()) * vertexCov(0, 1) +
168  px() * px() / (pt() * pt()) * vertexCov(1, 1));
169 }
size
Write out results.
uint8_t stopReason_
Stop Reason.
Definition: TrackBase.h:481
float chi2_
chi-squared
Definition: TrackBase.h:443
double t0() const
time at the reference point
Definition: TrackBase.h:726
void setQuality(const TrackQuality)
Definition: TrackBase.h:562
unsigned int index
index type
Definition: TrackBase.h:95
uint8_t quality_
track quality
Definition: TrackBase.h:475
uint8_t originalAlgorithm_
track algorithm
Definition: TrackBase.h:471
TrackBase()
default constructor
Definition: TrackBase.cc:69
TrackQuality
track quality
Definition: TrackBase.h:151
double dxyError() const
error on dxy
Definition: TrackBase.h:847
char charge_
electric charge
Definition: TrackBase.h:465
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:446
uint8_t stopReason() const
Definition: TrackBase.h:429
float covbetabeta_
Definition: TrackBase.h:440
ErrorD< N >::type type
Definition: Error.h:33
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:666
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:714
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:478
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:738
double beta() const
velocity at the reference point in natural units
Definition: TrackBase.h:732
static const std::string qualityNames[]
Definition: TrackBase.h:164
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:588
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:459
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:594
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:782
T sqrt(T t)
Definition: SSEVec.h:18
PerigeeCovarianceMatrix & fillCovariance(PerigeeCovarianceMatrix &v, const float *data)
double pt() const
track transverse momentum
Definition: TrackBase.h:660
double f[11][100]
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:134
float covt0t0_
errors for time and velocity (separate from cov for now)
Definition: TrackBase.h:440
virtual ~TrackBase()
virtual destructor
Definition: TrackBase.cc:124
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:28
float ndof_
number of degrees of freedom
Definition: TrackBase.h:462
static const std::string algoNames[]
Definition: TrackBase.h:148
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:549
double covt0t0() const
error on t0
Definition: TrackBase.h:871
float t0_
time at the reference point on track
Definition: TrackBase.h:449
fixed size matrix
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:437
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
CovarianceMatrix & fill(CovarianceMatrix &v) const
fill SMatrix
Definition: TrackBase.cc:129
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:146
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:452
int charge() const
track electric charge
Definition: TrackBase.h:606
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:468
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:80
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:672
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77