CMS 3D CMS Logo

L1TrackJetEmulatorProducer.cc
Go to the documentation of this file.
1 // Original Author: Samuel Edwin Leigh, Tyler Wu,
2 // Rutgers, the State University of New Jersey
3 //
4 // Rewritting/Improvements: George Karathanasis,
5 // georgios.karathanasis@cern.ch, CU Boulder
6 //
7 // Created: Wed, 01 Aug 2018 14:01:41 GMT
8 // Latest update: Nov 2022 (by GK)
9 //
10 // Track jets are clustered in a two-layer process, first by clustering in phi,
11 // then by clustering in eta. The code proceeds as following: putting all tracks// in a grid of eta vs phi space, and then cluster them. Finally we merge the cl
12 // usters when needed. The code is improved to use the same module between emula
13 // tion and simulation was also improved, with bug fixes and being faster.
14 // Introduction to object (p10-13):
15 // https://indico.cern.ch/event/791517/contributions/3341650/attachments/1818736/2973771/TrackBasedAlgos_L1TMadrid_MacDonald.pdf
16 // New and improved version: https://indico.cern.ch/event/1203796/contributions/5073056/attachments/2519806/4333006/trackjet_emu.pdf
17 
18 // L1T include files
26 
27 // system include files
39 
40 //own headers
41 #include "L1TrackJetClustering.h"
42 
43 //general
44 #include <ap_int.h>
45 
46 using namespace std;
47 using namespace edm;
48 using namespace l1t;
49 using namespace l1ttrackjet;
50 
52 public:
53  explicit L1TrackJetEmulatorProducer(const ParameterSet &);
54  ~L1TrackJetEmulatorProducer() override = default;
56  typedef vector<L1TTTrackType> L1TTTrackCollectionType;
58  static void fillDescriptions(ConfigurationDescriptions &descriptions);
59 
60 private:
61  void produce(Event &, const EventSetup &) override;
62 
63  // ----------member data ---------------------------
64 
65  std::vector<edm::Ptr<L1TTTrackType>> L1TrkPtrs_;
66  vector<int> tdtrk_;
67  const float trkZMax_;
68  const float trkPtMax_;
69  const float trkPtMin_;
70  const float trkEtaMax_;
71  const float nStubs4PromptChi2_;
72  const float nStubs5PromptChi2_;
73  const float nStubs4PromptBend_;
74  const float nStubs5PromptBend_;
75  const int trkNPSStubMin_;
77  const float lowpTJetThreshold_;
79  const float highpTJetThreshold_;
80  const int zBins_;
81  const int etaBins_;
82  const int phiBins_;
83  const double minTrkJetpT_;
84  const bool displaced_;
85  const float d0CutNStubs4_;
86  const float d0CutNStubs5_;
87  const float nStubs4DisplacedChi2_;
88  const float nStubs5DisplacedChi2_;
89  const float nStubs4DisplacedBend_;
90  const float nStubs5DisplacedBend_;
91  const int nDisplacedTracks_;
92  const float dzPVTrk_;
93 
94  float PVz;
95  float zStep_;
98 
100 
104 };
105 
106 //constructor
108  : trkZMax_(iConfig.getParameter<double>("trk_zMax")),
109  trkPtMax_(iConfig.getParameter<double>("trk_ptMax")),
110  trkPtMin_(iConfig.getParameter<double>("trk_ptMin")),
111  trkEtaMax_(iConfig.getParameter<double>("trk_etaMax")),
112  nStubs4PromptChi2_(iConfig.getParameter<double>("nStubs4PromptChi2")),
113  nStubs5PromptChi2_(iConfig.getParameter<double>("nStubs5PromptChi2")),
114  nStubs4PromptBend_(iConfig.getParameter<double>("nStubs4PromptBend")),
115  nStubs5PromptBend_(iConfig.getParameter<double>("nStubs5PromptBend")),
116  trkNPSStubMin_(iConfig.getParameter<int>("trk_nPSStubMin")),
117  lowpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
118  lowpTJetThreshold_(iConfig.getParameter<double>("lowpTJetThreshold")),
119  highpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
120  highpTJetThreshold_(iConfig.getParameter<double>("highpTJetThreshold")),
121  zBins_(iConfig.getParameter<int>("zBins")),
122  etaBins_(iConfig.getParameter<int>("etaBins")),
123  phiBins_(iConfig.getParameter<int>("phiBins")),
124  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
125  displaced_(iConfig.getParameter<bool>("displaced")),
126  d0CutNStubs4_(iConfig.getParameter<double>("d0_cutNStubs4")),
127  d0CutNStubs5_(iConfig.getParameter<double>("d0_cutNStubs5")),
128  nStubs4DisplacedChi2_(iConfig.getParameter<double>("nStubs4DisplacedChi2")),
129  nStubs5DisplacedChi2_(iConfig.getParameter<double>("nStubs5DisplacedChi2")),
130  nStubs4DisplacedBend_(iConfig.getParameter<double>("nStubs4DisplacedBend")),
131  nStubs5DisplacedBend_(iConfig.getParameter<double>("nStubs5DisplacedBend")),
132  nDisplacedTracks_(iConfig.getParameter<int>("nDisplacedTracks")),
133  dzPVTrk_(iConfig.getParameter<double>("MaxDzTrackPV")),
135  trackToken_(consumes<L1TTTrackRefCollectionType>(iConfig.getParameter<InputTag>("L1TrackInputTag"))),
136  PVtxToken_(consumes<l1t::VertexWordCollection>(iConfig.getParameter<InputTag>("L1PVertexInputTag"))) {
137  zStep_ = 2.0 * trkZMax_ / (zBins_ + 1); // added +1 in denom
138  etaStep_ = glbeta_intern(2.0 * trkEtaMax_ / etaBins_); //etaStep is the width of an etabin
139  phiStep_ = DoubleToBit(2.0 * (M_PI) / phiBins_,
140  TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit,
142 
143  if (displaced_)
144  produces<l1t::TkJetWordCollection>("L1TrackJetsExtended");
145  else
146  produces<l1t::TkJetWordCollection>("L1TrackJets");
147 }
148 
150  unique_ptr<l1t::TkJetWordCollection> L1TrackJetContainer(new l1t::TkJetWordCollection);
151  // Read inputs
152  const TrackerTopology &tTopo = iSetup.getData(tTopoToken_);
153 
155  iEvent.getByToken(trackToken_, TTTrackHandle);
156 
158  iEvent.getByToken(PVtxToken_, PVtx);
159  float PVz = (PVtx->at(0)).z0();
160 
161  L1TrkPtrs_.clear();
162  tdtrk_.clear();
163  // track selection
164  for (unsigned int this_l1track = 0; this_l1track < TTTrackHandle->size(); this_l1track++) {
165  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
166  float trk_pt = trkPtr->momentum().perp();
167  int trk_nstubs = (int)trkPtr->getStubRefs().size();
168  float trk_chi2dof = trkPtr->chi2Red();
169  float trk_bendchi2 = trkPtr->stubPtConsistency();
170  int trk_nPS = 0;
171  for (int istub = 0; istub < trk_nstubs; istub++) {
172  DetId detId(trkPtr->getStubRefs().at(istub)->getDetId());
173  if (detId.det() == DetId::Detector::Tracker) {
174  if ((detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3) ||
175  (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9))
176  trk_nPS++;
177  }
178  }
179  // selection tracks - supposed to happen on seperate module (kept for legacy/debug reasons)
180  if (trk_nPS < trkNPSStubMin_)
181  continue;
182  if (!TrackQualitySelection(trk_nstubs,
183  trk_chi2dof,
184  trk_bendchi2,
193  displaced_))
194  continue;
195  if (std::abs(PVz - trkPtr->z0()) > dzPVTrk_ && dzPVTrk_ > 0)
196  continue;
197  if (std::abs(trkPtr->z0()) > trkZMax_)
198  continue;
199  if (std::abs(trkPtr->momentum().eta()) > trkEtaMax_)
200  continue;
201  if (trk_pt < trkPtMin_)
202  continue;
203  L1TrkPtrs_.push_back(trkPtr);
204  }
205 
206  // if no tracks pass selection return empty containers
207  if (L1TrkPtrs_.empty()) {
208  if (displaced_)
209  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJetsExtended");
210  else
211  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJets");
212  return;
213  }
214 
216  mzb.znum = 0;
217  mzb.nclust = 0;
218  mzb.ht = 0;
219 
220  TrackJetEmulationEtaPhiBin epbins_default[phiBins_][etaBins_]; // create grid of phiBins
222  -1.0 * M_PI, TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit, TTTrack_TrackWord::stepPhi0);
223  for (int i = 0; i < phiBins_; ++i) {
225  for (int j = 0; j < etaBins_; ++j) {
226  epbins_default[i][j].phi = (phi + (phi + phiStep_)) / 2; // phi coord of bin center
227  epbins_default[i][j].eta = (eta + (eta + etaStep_)) / 2; // eta coord of bin center
228  eta += etaStep_;
229  } // for each etabin
230  phi += phiStep_;
231  } // for each phibin (finished creating epbins)
232 
233  // bins for z if multibin option is selected
234  std::vector<z0_intern> zmins, zmaxs;
235  for (int zbin = 0; zbin < zBins_; zbin++) {
236  zmins.push_back(DoubleToBit(
237  -1.0 * trkZMax_ + zStep_ * zbin, TTTrack_TrackWord::TrackBitWidths::kZ0Size, TTTrack_TrackWord::stepZ0));
238  zmaxs.push_back(DoubleToBit(-1.0 * trkZMax_ + zStep_ * zbin + 2 * zStep_,
239  TTTrack_TrackWord::TrackBitWidths::kZ0Size,
241  }
242 
243  // create vectors that hold clusters
244  std::vector<std::vector<TrackJetEmulationEtaPhiBin>> L1clusters;
245  L1clusters.reserve(phiBins_);
246  std::vector<TrackJetEmulationEtaPhiBin> L2clusters;
247 
248  // default (empty) grid
249  for (int i = 0; i < phiBins_; ++i) {
250  for (int j = 0; j < etaBins_; ++j) {
251  epbins_default[i][j].pTtot = 0;
252  epbins_default[i][j].used = false;
253  epbins_default[i][j].ntracks = 0;
254  epbins_default[i][j].nxtracks = 0;
255  epbins_default[i][j].trackidx.clear();
256  }
257  }
258 
259  // logic: loop over z bins find tracks in this bin and arrange them in a 2D eta-phi matrix
260  for (unsigned int zbin = 0; zbin < zmins.size(); ++zbin) {
261  // initialize matrices for every z bin
262  z0_intern zmin = zmins[zbin];
263  z0_intern zmax = zmaxs[zbin];
265  std::copy(&epbins_default[0][0], &epbins_default[0][0] + phiBins_ * etaBins_, &epbins[0][0]);
266 
267  //clear containers
268  L1clusters.clear();
269  L2clusters.clear();
270 
271  // fill grid
272  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
274  //-z0
275  z0_intern trkZ = L1TrkPtrs_[k]->getZ0Word();
276 
277  if (zmax < trkZ)
278  continue;
279  if (zmin > trkZ)
280  continue;
281  if (zbin == 0 && zmin == trkZ)
282  continue;
283 
284  //-pt
285  ap_uint<TTTrack_TrackWord::TrackBitWidths::kRinvSize - 1> ptEmulationBits = L1TrkPtrs_[k]->getRinvWord();
286  pt_intern trkpt;
287  trkpt.V = ptEmulationBits.range();
288 
289  //-eta
290  TTTrack_TrackWord::tanl_t etaEmulationBits = L1TrkPtrs_[k]->getTanlWord();
291  glbeta_intern trketa;
292  trketa.V = etaEmulationBits.range();
293 
294  //-phi
295  int Sector = L1TrkPtrs_[k]->phiSector();
296  double sector_phi_value = 0;
297  if (Sector < 5) {
298  sector_phi_value = 2.0 * M_PI * Sector / 9.0;
299  } else {
300  sector_phi_value = (-1.0 * M_PI + M_PI / 9.0 + (Sector - 5) * 2.0 * M_PI / 9.0);
301  }
302 
303  glbphi_intern trkphiSector = DoubleToBit(sector_phi_value,
304  TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit,
306  glbphi_intern local_phi = 0;
307  local_phi.V = L1TrkPtrs_[k]->getPhiWord();
308  glbphi_intern local_phi2 =
309  DoubleToBit(BitToDouble(local_phi, TTTrack_TrackWord::TrackBitWidths::kPhiSize, TTTrack_TrackWord::stepPhi0),
310  TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit,
312  glbphi_intern trkphi = local_phi2 + trkphiSector;
313 
314  //-d0
315  d0_intern abs_trkD0 = L1TrkPtrs_[k]->getD0Word();
316 
317  //-nstub
318  int trk_nstubs = (int)L1TrkPtrs_[k]->getStubRefs().size();
319 
320  // now fill the 2D grid with tracks
321  for (int i = 0; i < phiBins_; ++i) {
322  for (int j = 0; j < etaBins_; ++j) {
323  glbeta_intern eta_min = epbins[i][j].eta - etaStep_ / 2; //eta min
324  glbeta_intern eta_max = epbins[i][j].eta + etaStep_ / 2; //eta max
325  glbphi_intern phi_min = epbins[i][j].phi - phiStep_ / 2; //phi min
326  glbphi_intern phi_max = epbins[i][j].phi + phiStep_ / 2; //phi max
327  if ((trketa < eta_min) && j != 0)
328  continue;
329  if ((trketa > eta_max) && j != etaBins_ - 1)
330  continue;
331  if ((trkphi < phi_min) && i != 0)
332  continue;
333  if ((trkphi > phi_max) && i != phiBins_ - 1)
334  continue;
335 
336  if (trkpt < pt_intern(trkPtMax_))
337  epbins[i][j].pTtot += trkpt;
338  else
339  epbins[i][j].pTtot += pt_intern(trkPtMax_);
340  if ((abs_trkD0 >
341  DoubleToBit(d0CutNStubs5_, TTTrack_TrackWord::TrackBitWidths::kD0Size, TTTrack_TrackWord::stepD0) &&
342  trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
343  (abs_trkD0 >
344  DoubleToBit(d0CutNStubs4_, TTTrack_TrackWord::TrackBitWidths::kD0Size, TTTrack_TrackWord::stepD0) &&
345  trk_nstubs == 4 && d0CutNStubs4_ >= 0))
346  epbins[i][j].nxtracks += 1;
347 
348  epbins[i][j].trackidx.push_back(k);
349  ++epbins[i][j].ntracks;
350  } // for each etabin
351  } // for each phibin
352  } //end loop over tracks
353 
354  // first layer clustering - in eta using grid
355  for (int phibin = 0; phibin < phiBins_; ++phibin) {
356  L1clusters.push_back(L1_clustering<TrackJetEmulationEtaPhiBin, pt_intern, glbeta_intern, glbphi_intern>(
357  epbins[phibin], etaBins_, etaStep_));
358  }
359 
360  //second layer clustering in phi for all eta clusters
361  L2clusters = L2_clustering<TrackJetEmulationEtaPhiBin, pt_intern, glbeta_intern, glbphi_intern>(
362  L1clusters, phiBins_, phiStep_, etaStep_);
363 
364  // sum all cluster pTs in this zbin to find max
365  pt_intern sum_pt = 0;
366  for (unsigned int k = 0; k < L2clusters.size(); ++k) {
367  if (L2clusters[k].pTtot > pt_intern(highpTJetThreshold_) && L2clusters[k].ntracks < lowpTJetMinTrackMultiplicity_)
368  continue;
369  if (L2clusters[k].pTtot > pt_intern(highpTJetThreshold_) &&
370  L2clusters[k].ntracks < highpTJetMinTrackMultiplicity_)
371  continue;
372 
373  if (L2clusters[k].pTtot > pt_intern(minTrkJetpT_))
374  sum_pt += L2clusters[k].pTtot;
375  }
376  if (sum_pt < mzb.ht)
377  continue;
378  // if ht is larger than previous max, this is the new vertex zbin
379  mzb.ht = sum_pt;
380  mzb.znum = zbin;
381  mzb.clusters = L2clusters;
382  mzb.nclust = L2clusters.size();
383  mzb.zbincenter = (zmin + zmax) / 2.0;
384  } //zbin loop
385 
386  vector<edm::Ptr<L1TTTrackType>> L1TrackAssocJet;
387  for (unsigned int j = 0; j < mzb.clusters.size(); ++j) {
388  if (mzb.clusters[j].pTtot < pt_intern(trkPtMin_))
389  continue;
390 
392  TkJetWord::TkJetBitWidths::kGlbEtaSize,
393  TkJetWord::MAX_ETA / (1 << TkJetWord::TkJetBitWidths::kGlbEtaSize));
395  BitToDouble(mzb.clusters[j].phi, TTTrack_TrackWord::TrackBitWidths::kPhiSize + 4, TTTrack_TrackWord::stepPhi0),
396  TkJetWord::TkJetBitWidths::kGlbPhiSize,
397  (2. * std::abs(M_PI)) / (1 << TkJetWord::TkJetBitWidths::kGlbPhiSize));
398  l1t::TkJetWord::z0_t jetZ0 = 0;
399  l1t::TkJetWord::pt_t jetPt = mzb.clusters[j].pTtot;
400  l1t::TkJetWord::nt_t total_ntracks = mzb.clusters[j].ntracks;
401  l1t::TkJetWord::nx_t total_disptracks = mzb.clusters[j].nxtracks;
402  l1t::TkJetWord::dispflag_t dispflag = 0;
403  l1t::TkJetWord::tkjetunassigned_t unassigned = 0;
404 
405  if (total_disptracks > nDisplacedTracks_ || total_disptracks == nDisplacedTracks_)
406  dispflag = 1;
407  L1TrackAssocJet.clear();
408  for (unsigned int itrk = 0; itrk < mzb.clusters[j].trackidx.size(); itrk++)
409  L1TrackAssocJet.push_back(L1TrkPtrs_[mzb.clusters[j].trackidx[itrk]]);
410 
411  l1t::TkJetWord trkJet(jetPt, jetEta, jetPhi, jetZ0, total_ntracks, total_disptracks, dispflag, unassigned);
412 
413  L1TrackJetContainer->push_back(trkJet);
414  }
415 
416  std::sort(L1TrackJetContainer->begin(), L1TrackJetContainer->end(), [](auto &a, auto &b) { return a.pt() > b.pt(); });
417  if (displaced_)
418  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJetsExtended");
419  else
420  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJets");
421 }
422 
424  //The following says we do not know what parameters are allowed so do no validation
425  // Please change this to state exactly what you do use, even if it is no parameters
427  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("l1tTTTracksFromTrackletEmulation", "Level1TTTracks"));
428  desc.add<edm::InputTag>("L1PVertexInputTag", edm::InputTag("l1tVertexFinderEmulator", "L1VerticesEmulation"));
429  desc.add<double>("MaxDzTrackPV", 1.0);
430  desc.add<double>("trk_zMax", 15.0);
431  desc.add<double>("trk_ptMax", 200.0);
432  desc.add<double>("trk_ptMin", 3.0);
433  desc.add<double>("trk_etaMax", 2.4);
434  desc.add<double>("nStubs4PromptChi2", 5.0);
435  desc.add<double>("nStubs4PromptBend", 1.7);
436  desc.add<double>("nStubs5PromptChi2", 2.75);
437  desc.add<double>("nStubs5PromptBend", 3.5);
438  desc.add<int>("trk_nPSStubMin", -1);
439  desc.add<double>("minTrkJetpT", -1.0);
440  desc.add<int>("etaBins", 24);
441  desc.add<int>("phiBins", 27);
442  desc.add<int>("zBins", 1);
443  desc.add<double>("d0_cutNStubs4", -1);
444  desc.add<double>("d0_cutNStubs5", -1);
445  desc.add<int>("lowpTJetMinTrackMultiplicity", 2);
446  desc.add<double>("lowpTJetThreshold", 50.0);
447  desc.add<int>("highpTJetMinTrackMultiplicity", 3);
448  desc.add<double>("highpTJetThreshold", 100.0);
449  desc.add<bool>("displaced", false);
450  desc.add<double>("nStubs4DisplacedChi2", 5.0);
451  desc.add<double>("nStubs4DisplacedBend", 1.7);
452  desc.add<double>("nStubs5DisplacedChi2", 2.75);
453  desc.add<double>("nStubs5DisplacedBend", 3.5);
454  desc.add<int>("nDisplacedTracks", 2);
455  descriptions.add("l1tTrackJetsEmulator", desc);
456 }
457 
458 //define this as a plug-in
ap_uint< TTTrack_TrackWord::TrackBitWidths::kD0Size > d0_intern
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int tobLayer(const DetId &id) const
std::vector< l1t::TkJetWord > TkJetWordCollection
Definition: TkJetWord.h:191
T perp() const
Definition: PV3DBase.h:69
ap_int< kGlbEtaSize > glbeta_t
Definition: TkJetWord.h:59
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const unsigned int kExtraGlobalPhiBit
bool TrackQualitySelection(int trk_nstub, double trk_chi2, double trk_bendchi2, double nStubs4PromptBend_, double nStubs5PromptBend_, double nStubs4PromptChi2_, double nStubs5PromptChi2_, double nStubs4DisplacedBend_, double nStubs5DisplacedBend_, double nStubs4DisplacedChi2_, double nStubs5DisplacedChi2_, bool displaced_)
const unsigned int DoubleToBit(double value, unsigned int maxBits, double step)
T eta() const
Definition: PV3DBase.h:73
static constexpr double stepD0
ap_int< kGlbPhiSize > glbphi_t
Definition: TkJetWord.h:60
delete x;
Definition: CaloConfig.h:22
ap_ufixed< kPtSize, kPtMagSize, AP_TRN, AP_SAT > pt_t
Definition: TkJetWord.h:58
ap_uint< kDispFlagSize > dispflag_t
Definition: TkJetWord.h:64
ap_uint< TkJetBitWidths::kUnassignedSize > tkjetunassigned_t
Definition: TkJetWord.h:65
ap_int< kZ0Size > z0_t
Definition: TkJetWord.h:61
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
const double BitToDouble(unsigned int bits, unsigned int maxBits, double step)
ap_fixed< TTTrack_TrackWord::TrackBitWidths::kTanlSize, ETA_INTPART_BITS, AP_TRN, AP_SAT > glbeta_intern
edm::RefVector< L1TTTrackCollectionType > L1TTTrackRefCollectionType
int iEvent
Definition: GenABIO.cc:224
L1TrackJetEmulatorProducer(const ParameterSet &)
ap_ufixed< TTTrack_TrackWord::TrackBitWidths::kRinvSize - 1, PT_INTPART_BITS, AP_TRN, AP_SAT > pt_intern
vector< L1TTTrackType > L1TTTrackCollectionType
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(Event &, const EventSetup &) override
ap_int< TTTrack_TrackWord::TrackBitWidths::kZ0Size > z0_intern
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< VertexWord > VertexWordCollection
Definition: VertexWord.h:197
static constexpr auto TOB
std::vector< TrackJetEmulationEtaPhiBin > clusters
ap_uint< kXtSize > nx_t
Definition: TkJetWord.h:63
std::vector< edm::Ref< edmNew::DetSetVector< TTStub< T > >, TTStub< T > > > getStubRefs() const
Track components.
Definition: TTTrack.h:93
const EDGetTokenT< L1TTTrackRefCollectionType > trackToken_
#define M_PI
Definition: DetId.h:17
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ap_int< TTTrack_TrackWord::TrackBitWidths::kPhiSize+kExtraGlobalPhiBit > glbphi_intern
std::vector< edm::Ptr< L1TTTrackType > > L1TrkPtrs_
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
double stubPtConsistency() const
StubPtConsistency.
Definition: TTTrack.h:417
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(ConfigurationDescriptions &descriptions)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
static constexpr double stepPhi0
const EDGetTokenT< l1t::VertexWordCollection > PVtxToken_
HLT enums.
double chi2Red() const
Chi2 reduced.
Definition: TTTrack.h:359
ap_uint< kNtSize > nt_t
Definition: TkJetWord.h:62
ap_uint< TrackBitWidths::kTanlSize > tanl_t
double a
Definition: hdecay.h:121
double z0() const
Track z0.
Definition: TTTrack.h:330
unsigned int tidRing(const DetId &id) const
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
static constexpr double stepZ0
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511