CMS 3D CMS Logo

L1track3D.h
Go to the documentation of this file.
1 #ifndef L1Trigger_TrackFindingTMTT_L1track3D_h
2 #define L1Trigger_TrackFindingTMTT_L1track3D_h
3 
12 
13 #include <vector>
14 #include <string>
15 #include <unordered_set>
16 #include <utility>
17 
18 //=== L1 track candidate found in 3 dimensions.
19 //=== Gives access to all stubs on track and to its 3D helix parameters.
20 //=== Also calculates & gives access to associated truth particle (Tracking Particle) if any.
21 
22 namespace tmtt {
23 
24  class L1track3D : public L1trackBase {
25  public:
26  // Seeding layers of tracklet pattern reco.
28 
29  public:
30  L1track3D(const Settings* settings,
31  const std::vector<Stub*>& stubs,
32  std::pair<unsigned int, unsigned int> cellLocationHT,
33  std::pair<float, float> helixRphi,
34  std::pair<float, float> helixRz,
35  float helixD0,
36  unsigned int iPhiSec,
37  unsigned int iEtaReg,
38  unsigned int optoLinkID,
39  bool mergedHTcell)
40  : L1trackBase(),
41  settings_(settings),
42  stubs_(stubs),
43  stubsConst_(stubs_.begin(), stubs_.end()),
47  helixD0_(helixD0),
53  seedPS_(999) {
54  nLayers_ = Utility::countLayers(settings, stubs_); // Count tracker layers these stubs are in
56  stubs_,
58  matchedStubs_); // Find associated truth particle & calculate info about match.
59  }
60 
61  // TMTT tracking: constructor
62 
63  L1track3D(const Settings* settings,
64  const std::vector<Stub*>& stubs,
65  std::pair<unsigned int, unsigned int> cellLocationHT,
66  std::pair<float, float> helixRphi,
67  std::pair<float, float> helixRz,
68  unsigned int iPhiSec,
69  unsigned int iEtaReg,
70  unsigned int optoLinkID,
71  bool mergedHTcell)
72  : L1track3D(
74 
75  ~L1track3D() override = default;
76 
77  //--- Set/get optional info for tracklet tracks.
78 
79  // Tracklet seeding layer pair (from Tracklet::calcSeedIndex())
82 
83  // Tracklet seed stub pair uses PS modules (from FPGATracket::PSseed())
84  void setSeedPS(unsigned int seedPS) { seedPS_ = seedPS; }
85  unsigned int seedPS() const { return seedPS_; }
86 
87  // Best stub (stub with smallest Phi residual in each layer/disk)
88  void setBestStubs(std::unordered_set<const Stub*> bestStubs) { bestStubs_ = bestStubs; }
89  std::unordered_set<const Stub*> bestStubs() const { return bestStubs_; }
90 
91  //--- Get information about the reconstructed track.
92 
93  // Get stubs on track candidate.
94  const std::vector<const Stub*>& stubsConst() const override { return stubsConst_; }
95  const std::vector<Stub*>& stubs() const override { return stubs_; }
96  // Get number of stubs on track candidate.
97  unsigned int numStubs() const override { return stubs_.size(); }
98  // Get number of tracker layers these stubs are in.
99  unsigned int numLayers() const override { return nLayers_; }
100  // Get cell location of track candidate in r-phi Hough Transform array in units of bin number.
101  std::pair<unsigned int, unsigned int> cellLocationHT() const override { return cellLocationHT_; }
102  // The two conventionally agreed track helix parameters relevant in r-phi plane. i.e. (q/Pt, phi0)
103  std::pair<float, float> helixRphi() const { return helixRphi_; }
104  // The two conventionally agreed track helix parameters relevant in r-z plane. i.e. (z0, tan_lambda)
105  std::pair<float, float> helixRz() const { return helixRz_; }
106 
107  //--- Return chi variables, (both digitized & undigitized), which are the stub coords. relative to track.
108 
109  std::vector<float> chiPhi() {
110  std::vector<float> result;
111  for (const Stub* s : stubs_) {
112  float chi_phi = reco::deltaPhi(s->phi(), this->phi0() - s->r() * this->qOverPt() * settings_->invPtToDphi());
113  result.push_back(chi_phi);
114  }
115  return result;
116  }
117 
118  std::vector<int> chiPhiDigi() {
119  std::vector<int> result;
120  const float phiMult = pow(2, settings_->phiSBits()) / settings_->phiSRange();
121  for (const float& chi_phi : this->chiPhi()) {
122  int iDigi_chi_phi = floor(chi_phi * phiMult);
123  result.push_back(iDigi_chi_phi);
124  }
125  return result;
126  }
127 
128  std::vector<float> chiZ() {
129  std::vector<float> result;
130  for (const Stub* s : stubs_) {
131  float chi_z = s->z() - (this->z0() + s->r() * this->tanLambda());
132  result.push_back(chi_z);
133  }
134  return result;
135  }
136 
137  std::vector<int> chiZDigi() {
138  std::vector<int> result;
139  const float zMult = pow(2, settings_->zBits()) / settings_->zRange();
140  for (const float& chi_z : this->chiZ()) {
141  int iDigi_chi_z = floor(chi_z * zMult);
142  result.push_back(iDigi_chi_z);
143  }
144  return result;
145  }
146 
147  //--- User-friendlier access to the helix parameters.
148 
149  float qOverPt() const override { return helixRphi_.first; }
150  float charge() const { return (this->qOverPt() > 0 ? 1 : -1); }
151  float invPt() const { return std::abs(this->qOverPt()); }
152  // Protect pt against 1/pt = 0.
153  float pt() const {
154  constexpr float small = 1.0e-6;
155  return 1. / (small + this->invPt());
156  }
157  float d0() const { return helixD0_; } // Hough transform assumes d0 = 0.
158  float phi0() const override { return helixRphi_.second; }
159  float z0() const { return helixRz_.first; }
160  float tanLambda() const { return helixRz_.second; }
161  float theta() const { return atan2(1., this->tanLambda()); } // Use atan2 to ensure 0 < theta < pi.
162  float eta() const { return -log(tan(0.5 * this->theta())); }
163 
164  // Phi and z coordinates at which track crosses "chosenR" values used by r-phi HT and rapidity sectors respectively.
165  float phiAtChosenR() const {
166  return reco::deltaPhi(this->phi0() - (settings_->invPtToDphi() * settings_->chosenRofPhi()) * this->qOverPt(),
167  0.);
168  }
169  float zAtChosenR() const {
170  return (this->z0() + (settings_->chosenRofZ()) * this->tanLambda());
171  } // neglects transverse impact parameter & track curvature.
172 
173  //--- Get phi sector and eta region used by track finding code that this track is in.
174  unsigned int iPhiSec() const override { return iPhiSec_; }
175  unsigned int iEtaReg() const override { return iEtaReg_; }
176 
177  //--- Opto-link ID used to send this track from HT to Track Fitter
178  unsigned int optoLinkID() const override { return optoLinkID_; }
179 
180  //--- Was this track produced from a marged HT cell (e.g. 2x2)?
181  bool mergedHTcell() const { return mergedHTcell_; }
182 
183  //--- Get information about its association (if any) to a truth Tracking Particle.
184 
185  // Get best matching tracking particle (=nullptr if none).
186  const TP* matchedTP() const override { return matchedTP_; }
187  // Get the matched stubs with this Tracking Particle
188  const std::vector<const Stub*>& matchedStubs() const override { return matchedStubs_; }
189  // Get number of matched stubs with this Tracking Particle
190  unsigned int numMatchedStubs() const override { return matchedStubs_.size(); }
191  // Get number of tracker layers with matched stubs with this Tracking Particle
192  unsigned int numMatchedLayers() const override { return nMatchedLayers_; }
193  // Get purity of stubs on track candidate (i.e. fraction matching best Tracking Particle)
194  float purity() const { return numMatchedStubs() / float(numStubs()); }
195 
196  //--- For debugging purposes.
197 
198  // Remove incorrect stubs from the track using truth information.
199  // Also veto tracks where the HT cell estimated from the true helix parameters is inconsistent with the cell the HT found the track in, (since probable duplicates).
200  // Also veto tracks that match a truth particle not used for the algo efficiency measurement.
201  // Return a boolean indicating if the track should be kept. (i.e. Is genuine & non-duplicate).
202  bool cheat() {
203  bool keep = false;
204 
205  std::vector<Stub*> stubsSel;
206  if (matchedTP_ != nullptr) { // Genuine track
207  for (Stub* s : stubs_) {
208  const TP* tp = s->assocTP();
209  if (tp != nullptr) {
210  if (matchedTP_->index() == tp->index()) {
211  stubsSel.push_back(s); // This stub was produced by same truth particle as rest of track, so keep it.
212  }
213  }
214  }
215  }
216  stubs_ = stubsSel;
217  stubsConst_ = std::vector<const Stub*>(stubs_.begin(), stubs_.end());
218 
219  nLayers_ = Utility::countLayers(settings_, stubs_); // Count tracker layers these stubs are in
221  stubs_,
223  matchedStubs_); // Find associated truth particle & calculate info about match.
224 
225  bool genuine = (matchedTP_ != nullptr);
226 
227  if (genuine && matchedTP_->useForAlgEff()) {
228  Sector secTmp(settings_, iPhiSec_, iEtaReg_);
229  HTrphi htRphiTmp(settings_, iPhiSec_, iEtaReg_, secTmp.etaMin(), secTmp.etaMax(), secTmp.phiCentre());
230  std::pair<unsigned int, unsigned int> trueCell = htRphiTmp.trueCell(matchedTP_);
231 
232  std::pair<unsigned int, unsigned int> htCell = this->cellLocationHT();
233  bool consistent = (htCell == trueCell); // If true, track is probably not a duplicate.
234  if (mergedHTcell_) {
235  // If this is a merged cell, check other elements of merged cell.
236  std::pair<unsigned int, unsigned int> htCell10(htCell.first + 1, htCell.second);
237  std::pair<unsigned int, unsigned int> htCell01(htCell.first, htCell.second + 1);
238  std::pair<unsigned int, unsigned int> htCell11(htCell.first + 1, htCell.second + 1);
239  if (htCell10 == trueCell)
240  consistent = true;
241  if (htCell01 == trueCell)
242  consistent = true;
243  if (htCell11 == trueCell)
244  consistent = true;
245  }
246  if (consistent)
247  keep = true;
248  }
249 
250  return keep; // Indicate if track should be kept.
251  }
252 
253  private:
254  //--- Configuration parameters
256 
257  //--- Information about the reconstructed track.
258  std::vector<Stub*> stubs_;
259  std::vector<const Stub*> stubsConst_;
260  std::unordered_set<const Stub*> bestStubs_;
261  unsigned int nLayers_;
262  std::pair<unsigned int, unsigned int> cellLocationHT_;
263  std::pair<float, float> helixRphi_;
264  std::pair<float, float> helixRz_;
265  float helixD0_;
266  unsigned int iPhiSec_;
267  unsigned int iEtaReg_;
268  unsigned int optoLinkID_;
270 
271  //--- Optional info used for tracklet tracks.
273  unsigned int seedPS_;
274 
275  //--- Information about its association (if any) to a truth Tracking Particle.
276  const TP* matchedTP_;
277  std::vector<const Stub*> matchedStubs_;
278  unsigned int nMatchedLayers_;
279  };
280 
281 } // namespace tmtt
282 
283 #endif
std::pair< unsigned int, unsigned int > cellLocationHT_
Definition: L1track3D.h:262
unsigned int iEtaReg() const override
Definition: L1track3D.h:175
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
TrackletSeedType seedLayerType() const
Definition: L1track3D.h:81
float phi0() const override
Definition: L1track3D.h:158
double invPtToDphi() const
Definition: Settings.h:397
unsigned int numMatchedStubs() const override
Definition: L1track3D.h:190
float qOverPt() const override
Definition: L1track3D.h:149
void setSeedLayerType(unsigned int seedLayerType)
Definition: L1track3D.h:80
std::vector< int > chiZDigi()
Definition: L1track3D.h:137
bool useForAlgEff() const
Definition: TP.h:89
TrackletSeedType seedLayerType_
Definition: L1track3D.h:272
bool mergedHTcell_
Definition: L1track3D.h:269
unsigned int optoLinkID_
Definition: L1track3D.h:268
std::pair< float, float > helixRz() const
Definition: L1track3D.h:105
const Settings * settings_
Definition: L1track3D.h:255
unsigned int nMatchedLayers_
Definition: L1track3D.h:278
double phiSRange() const
Definition: Settings.h:84
const TP * matchedTP() const override
Definition: L1track3D.h:186
float invPt() const
Definition: L1track3D.h:151
float z0() const
Definition: L1track3D.h:159
unsigned int seedPS() const
Definition: L1track3D.h:85
float purity() const
Definition: L1track3D.h:194
const std::vector< const Stub * > & stubsConst() const override
Definition: L1track3D.h:94
unsigned int numMatchedLayers() const override
Definition: L1track3D.h:192
const std::vector< Stub * > & stubs() const override
Definition: L1track3D.h:95
float pt() const
Definition: L1track3D.h:153
std::unordered_set< const Stub * > bestStubs_
Definition: L1track3D.h:260
std::pair< unsigned int, unsigned int > cellLocationHT() const override
Definition: L1track3D.h:101
float etaMax() const
Definition: Sector.h:31
float tanLambda() const
Definition: L1track3D.h:160
double zRange() const
Definition: Settings.h:88
std::pair< float, float > helixRphi() const
Definition: L1track3D.h:103
double chosenRofZ() const
Definition: Settings.h:127
void setBestStubs(std::unordered_set< const Stub *> bestStubs)
Definition: L1track3D.h:88
~L1track3D() override=default
bool mergedHTcell() const
Definition: L1track3D.h:181
unsigned int phiSBits() const
Definition: Settings.h:83
const std::vector< const Stub * > & matchedStubs() const override
Definition: L1track3D.h:188
unsigned int iPhiSec() const override
Definition: L1track3D.h:174
Definition: TP.h:23
std::pair< unsigned int, unsigned int > trueCell(const TP *tp) const override
Definition: HTrphi.cc:473
unsigned int index() const
Definition: TP.h:39
L1track3D(const Settings *settings, const std::vector< Stub *> &stubs, std::pair< unsigned int, unsigned int > cellLocationHT, std::pair< float, float > helixRphi, std::pair< float, float > helixRz, unsigned int iPhiSec, unsigned int iEtaReg, unsigned int optoLinkID, bool mergedHTcell)
Definition: L1track3D.h:63
unsigned int numLayers() const override
Definition: L1track3D.h:99
unsigned int zBits() const
Definition: Settings.h:87
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
const TP * matchingTP(const Settings *settings, const std::vector< const Stub *> &vstubs, unsigned int &nMatchedLayersBest, std::vector< const Stub *> &matchedStubsBest)
Definition: Utility.cc:63
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int iPhiSec_
Definition: L1track3D.h:266
std::vector< int > chiPhiDigi()
Definition: L1track3D.h:118
std::vector< float > chiPhi()
Definition: L1track3D.h:109
double chosenRofPhi() const
Definition: Settings.h:112
std::unordered_set< const Stub * > bestStubs() const
Definition: L1track3D.h:89
const TP * matchedTP_
Definition: L1track3D.h:276
float phiCentre() const
Definition: Sector.h:29
void setSeedPS(unsigned int seedPS)
Definition: L1track3D.h:84
std::vector< Stub * > stubs_
Definition: L1track3D.h:258
=== This is the base class for the linearised chi-squared track fit algorithms.
Definition: Array2D.h:16
float theta() const
Definition: L1track3D.h:161
std::vector< const Stub * > matchedStubs_
Definition: L1track3D.h:277
float etaMin() const
Definition: Sector.h:30
float eta() const
Definition: L1track3D.h:162
float zAtChosenR() const
Definition: L1track3D.h:169
unsigned int iEtaReg_
Definition: L1track3D.h:267
unsigned int optoLinkID() const override
Definition: L1track3D.h:178
unsigned int numStubs() const override
Definition: L1track3D.h:97
std::pair< float, float > helixRphi_
Definition: L1track3D.h:263
float phiAtChosenR() const
Definition: L1track3D.h:165
std::vector< const Stub * > stubsConst_
Definition: L1track3D.h:259
float charge() const
Definition: L1track3D.h:150
unsigned int nLayers_
Definition: L1track3D.h:261
float d0() const
Definition: L1track3D.h:157
unsigned int countLayers(const Settings *settings, const std::vector< const Stub *> &stubs, bool disableReducedLayerID=false, bool onlyPS=false)
Definition: Utility.cc:25
unsigned int seedPS_
Definition: L1track3D.h:273
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
std::vector< float > chiZ()
Definition: L1track3D.h:128
std::pair< float, float > helixRz_
Definition: L1track3D.h:264
L1track3D(const Settings *settings, const std::vector< Stub *> &stubs, std::pair< unsigned int, unsigned int > cellLocationHT, std::pair< float, float > helixRphi, std::pair< float, float > helixRz, float helixD0, unsigned int iPhiSec, unsigned int iEtaReg, unsigned int optoLinkID, bool mergedHTcell)
Definition: L1track3D.h:30