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 // Claire Savard (claire.savard@colorado.edu)
7 //
8 // Created: Wed, 01 Aug 2018 14:01:41 GMT
9 // Latest update: Nov 2023 (by CS)
10 //
11 // Track jets are clustered in a two-layer process, first by clustering in phi,
12 // 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
13 // usters when needed. The code is improved to use the same module between emula
14 // tion and simulation was also improved, with bug fixes and being faster.
15 // Introduction to object (p10-13):
16 // https://indico.cern.ch/event/791517/contributions/3341650/attachments/1818736/2973771/TrackBasedAlgos_L1TMadrid_MacDonald.pdf
17 // New and improved version: https://indico.cern.ch/event/1203796/contributions/5073056/attachments/2519806/4333006/trackjet_emu.pdf
18 
19 // L1T include files
27 
28 // system include files
40 
41 //own headers
42 #include "L1TrackJetClustering.h"
43 
44 //general
45 #include <ap_int.h>
46 
47 using namespace std;
48 using namespace edm;
49 using namespace l1t;
50 using namespace l1ttrackjet;
51 
53 public:
54  explicit L1TrackJetEmulatorProducer(const ParameterSet &);
55  ~L1TrackJetEmulatorProducer() override = default;
57  typedef vector<L1TTTrackType> L1TTTrackCollectionType;
59  static void fillDescriptions(ConfigurationDescriptions &descriptions);
60 
61 private:
62  void produce(Event &, const EventSetup &) override;
63 
64  // ----------member data ---------------------------
65 
66  std::vector<edm::Ptr<L1TTTrackType>> L1TrkPtrs_;
67  const float trkZMax_;
68  const float trkPtMax_;
69  const float trkEtaMax_;
71  const float lowpTJetThreshold_;
73  const float highpTJetThreshold_;
74  const int zBins_;
75  const int etaBins_;
76  const int phiBins_;
77  const double minTrkJetpT_;
78  const bool displaced_;
79  const float d0CutNStubs4_;
80  const float d0CutNStubs5_;
81  const int nDisplacedTracks_;
82 
83  float zStep_;
86 
88 
90 };
91 
92 //constructor
94  : trkZMax_(iConfig.getParameter<double>("trk_zMax")),
95  trkPtMax_(iConfig.getParameter<double>("trk_ptMax")),
96  trkEtaMax_(iConfig.getParameter<double>("trk_etaMax")),
97  lowpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
98  lowpTJetThreshold_(iConfig.getParameter<double>("lowpTJetThreshold")),
99  highpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
100  highpTJetThreshold_(iConfig.getParameter<double>("highpTJetThreshold")),
101  zBins_(iConfig.getParameter<int>("zBins")),
102  etaBins_(iConfig.getParameter<int>("etaBins")),
103  phiBins_(iConfig.getParameter<int>("phiBins")),
104  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
105  displaced_(iConfig.getParameter<bool>("displaced")),
106  d0CutNStubs4_(iConfig.getParameter<double>("d0_cutNStubs4")),
107  d0CutNStubs5_(iConfig.getParameter<double>("d0_cutNStubs5")),
108  nDisplacedTracks_(iConfig.getParameter<int>("nDisplacedTracks")),
109  trackToken_(consumes<L1TTTrackRefCollectionType>(iConfig.getParameter<InputTag>("L1TrackInputTag"))) {
110  zStep_ = 2.0 * trkZMax_ / (zBins_ + 1); // added +1 in denom
111  etaStep_ = glbeta_intern(2.0 * trkEtaMax_ / etaBins_); //etaStep is the width of an etabin
112  phiStep_ = DoubleToBit(2.0 * (M_PI) / phiBins_,
113  TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit,
115 
116  if (displaced_)
117  produces<l1t::TkJetWordCollection>("L1TrackJetsExtended");
118  else
119  produces<l1t::TkJetWordCollection>("L1TrackJets");
120 }
121 
123  unique_ptr<l1t::TkJetWordCollection> L1TrackJetContainer(new l1t::TkJetWordCollection);
124 
125  // L1 tracks
127  iEvent.getByToken(trackToken_, TTTrackHandle);
128 
129  L1TrkPtrs_.clear();
130  // track selection
131  for (unsigned int this_l1track = 0; this_l1track < TTTrackHandle->size(); this_l1track++) {
132  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
133  L1TrkPtrs_.push_back(trkPtr);
134  }
135 
136  // if no tracks pass selection return empty containers
137  if (L1TrkPtrs_.empty()) {
138  if (displaced_)
139  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJetsExtended");
140  else
141  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJets");
142  return;
143  }
144 
146  mzb.znum = 0;
147  mzb.nclust = 0;
148  mzb.ht = 0;
149 
150  TrackJetEmulationEtaPhiBin epbins_default[phiBins_][etaBins_]; // create grid of phiBins
152  -1.0 * M_PI, TTTrack_TrackWord::TrackBitWidths::kPhiSize + kExtraGlobalPhiBit, TTTrack_TrackWord::stepPhi0);
153  for (int i = 0; i < phiBins_; ++i) {
155  for (int j = 0; j < etaBins_; ++j) {
156  epbins_default[i][j].phi = (phi + (phi + phiStep_)) / 2; // phi coord of bin center
157  epbins_default[i][j].eta = (eta + (eta + etaStep_)) / 2; // eta coord of bin center
158  eta += etaStep_;
159  } // for each etabin
160  phi += phiStep_;
161  } // for each phibin (finished creating epbins)
162 
163  // bins for z if multibin option is selected
164  std::vector<z0_intern> zmins, zmaxs;
165  for (int zbin = 0; zbin < zBins_; zbin++) {
166  zmins.push_back(DoubleToBit(
167  -1.0 * trkZMax_ + zStep_ * zbin, TTTrack_TrackWord::TrackBitWidths::kZ0Size, TTTrack_TrackWord::stepZ0));
168  zmaxs.push_back(DoubleToBit(-1.0 * trkZMax_ + zStep_ * zbin + 2 * zStep_,
169  TTTrack_TrackWord::TrackBitWidths::kZ0Size,
171  }
172 
173  // create vectors that hold clusters
174  std::vector<std::vector<TrackJetEmulationEtaPhiBin>> L1clusters;
175  L1clusters.reserve(phiBins_);
176  std::vector<TrackJetEmulationEtaPhiBin> L2clusters;
177 
178  // default (empty) grid
179  for (int i = 0; i < phiBins_; ++i) {
180  for (int j = 0; j < etaBins_; ++j) {
181  epbins_default[i][j].pTtot = 0;
182  epbins_default[i][j].used = false;
183  epbins_default[i][j].ntracks = 0;
184  epbins_default[i][j].nxtracks = 0;
185  epbins_default[i][j].trackidx.clear();
186  }
187  }
188 
189  //Begin Firmware-style clustering
190  // logic: loop over z bins find tracks in this bin and arrange them in a 2D eta-phi matrix
191  for (unsigned int zbin = 0; zbin < zmins.size(); ++zbin) {
192  // initialize matrices for every z bin
193  z0_intern zmin = zmins[zbin];
194  z0_intern zmax = zmaxs[zbin];
195 
197 
198  std::copy(&epbins_default[0][0], &epbins_default[0][0] + phiBins_ * etaBins_, &epbins[0][0]);
199 
200  L1clusters.clear();
201  L2clusters.clear();
202  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
203  z0_intern trkZ = L1TrkPtrs_[k]->getZ0Word();
204 
205  if (zmax < trkZ)
206  continue;
207  if (zmin > trkZ)
208  continue;
209  if (zbin == 0 && zmin == trkZ)
210  continue;
211 
212  // Pt
213  ap_uint<TTTrack_TrackWord::TrackBitWidths::kRinvSize - 1> ptEmulationBits = L1TrkPtrs_[k]->getRinvWord();
214  pt_intern trkpt;
215  trkpt.V = ptEmulationBits.range();
216 
217  // d0
218  d0_intern abs_trkD0 = L1TrkPtrs_[k]->getD0Word();
219 
220  // nstubs
221  int trk_nstubs = (int)L1TrkPtrs_[k]->getStubRefs().size();
222 
223  // Phi bin
224  int i = phi_bin_firmwareStyle(L1TrkPtrs_[k]->phiSector(),
225  L1TrkPtrs_[k]->getPhiWord()); //Function defined in L1TrackJetClustering.h
226 
227  // Eta bin
228  int j = eta_bin_firmwareStyle(L1TrkPtrs_[k]->getTanlWord()); //Function defined in L1TrackJetClustering.h
229 
230  if (trkpt < pt_intern(trkPtMax_))
231  epbins[i][j].pTtot += trkpt;
232  else
233  epbins[i][j].pTtot += pt_intern(trkPtMax_);
234  if ((abs_trkD0 >
235  DoubleToBit(d0CutNStubs5_, TTTrack_TrackWord::TrackBitWidths::kD0Size, TTTrack_TrackWord::stepD0) &&
236  trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
237  (abs_trkD0 >
238  DoubleToBit(d0CutNStubs4_, TTTrack_TrackWord::TrackBitWidths::kD0Size, TTTrack_TrackWord::stepD0) &&
239  trk_nstubs == 4 && d0CutNStubs4_ >= 0))
240  epbins[i][j].nxtracks += 1;
241 
242  epbins[i][j].trackidx.push_back(k);
243  ++epbins[i][j].ntracks;
244  }
245  //End Firmware style clustering
246 
247  // first layer clustering - in eta using grid
248  for (int phibin = 0; phibin < phiBins_; ++phibin) {
249  L1clusters.push_back(L1_clustering<TrackJetEmulationEtaPhiBin, pt_intern, glbeta_intern, glbphi_intern>(
250  epbins[phibin], etaBins_, etaStep_));
251  }
252 
253  //second layer clustering in phi for all eta clusters
254  L2clusters = L2_clustering<TrackJetEmulationEtaPhiBin, pt_intern, glbeta_intern, glbphi_intern>(
255  L1clusters, phiBins_, phiStep_, etaStep_);
256 
257  // sum all cluster pTs in this zbin to find max
258  pt_intern sum_pt = 0;
259  for (unsigned int k = 0; k < L2clusters.size(); ++k) {
260  if (L2clusters[k].pTtot > pt_intern(highpTJetThreshold_) && L2clusters[k].ntracks < lowpTJetMinTrackMultiplicity_)
261  continue;
262  if (L2clusters[k].pTtot > pt_intern(highpTJetThreshold_) &&
263  L2clusters[k].ntracks < highpTJetMinTrackMultiplicity_)
264  continue;
265 
266  if (L2clusters[k].pTtot > pt_intern(minTrkJetpT_))
267  sum_pt += L2clusters[k].pTtot;
268  }
269  if (sum_pt < mzb.ht)
270  continue;
271  // if ht is larger than previous max, this is the new vertex zbin
272  mzb.ht = sum_pt;
273  mzb.znum = zbin;
274  mzb.clusters = L2clusters;
275  mzb.nclust = L2clusters.size();
276  mzb.zbincenter = (zmin + zmax) / 2.0;
277  } //zbin loop
278 
279  vector<edm::Ptr<L1TTTrackType>> L1TrackAssocJet;
280  for (unsigned int j = 0; j < mzb.clusters.size(); ++j) {
282  TkJetWord::TkJetBitWidths::kGlbEtaSize,
283  TkJetWord::MAX_ETA / (1 << TkJetWord::TkJetBitWidths::kGlbEtaSize));
285  BitToDouble(mzb.clusters[j].phi, TTTrack_TrackWord::TrackBitWidths::kPhiSize + 4, TTTrack_TrackWord::stepPhi0),
286  TkJetWord::TkJetBitWidths::kGlbPhiSize,
287  (2. * std::abs(M_PI)) / (1 << TkJetWord::TkJetBitWidths::kGlbPhiSize));
288  l1t::TkJetWord::z0_t jetZ0 = 0;
289  l1t::TkJetWord::pt_t jetPt = mzb.clusters[j].pTtot;
290  l1t::TkJetWord::nt_t total_ntracks = mzb.clusters[j].ntracks;
291  l1t::TkJetWord::nx_t total_disptracks = mzb.clusters[j].nxtracks;
292  l1t::TkJetWord::dispflag_t dispflag = 0;
293  l1t::TkJetWord::tkjetunassigned_t unassigned = 0;
294 
295  if (total_disptracks >= nDisplacedTracks_)
296  dispflag = 1;
297  L1TrackAssocJet.clear();
298  for (unsigned int itrk = 0; itrk < mzb.clusters[j].trackidx.size(); itrk++)
299  L1TrackAssocJet.push_back(L1TrkPtrs_[mzb.clusters[j].trackidx[itrk]]);
300 
301  l1t::TkJetWord trkJet(jetPt, jetEta, jetPhi, jetZ0, total_ntracks, total_disptracks, dispflag, unassigned);
302 
303  L1TrackJetContainer->push_back(trkJet);
304  }
305 
306  std::sort(L1TrackJetContainer->begin(), L1TrackJetContainer->end(), [](auto &a, auto &b) { return a.pt() > b.pt(); });
307  if (displaced_)
308  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJetsExtended");
309  else
310  iEvent.put(std::move(L1TrackJetContainer), "L1TrackJets");
311 }
312 
314  //The following says we do not know what parameters are allowed so do no validation
315  // Please change this to state exactly what you do use, even if it is no parameters
317  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("l1tTTTracksFromTrackletEmulation", "Level1TTTracks"));
318  desc.add<double>("trk_zMax", 15.0);
319  desc.add<double>("trk_ptMax", 200.0);
320  desc.add<double>("trk_etaMax", 2.4);
321  desc.add<double>("minTrkJetpT", -1.0);
322  desc.add<int>("etaBins", 24);
323  desc.add<int>("phiBins", 27);
324  desc.add<int>("zBins", 1);
325  desc.add<double>("d0_cutNStubs4", -1);
326  desc.add<double>("d0_cutNStubs5", -1);
327  desc.add<int>("lowpTJetMinTrackMultiplicity", 2);
328  desc.add<double>("lowpTJetThreshold", 50.0);
329  desc.add<int>("highpTJetMinTrackMultiplicity", 3);
330  desc.add<double>("highpTJetThreshold", 100.0);
331  desc.add<bool>("displaced", false);
332  desc.add<int>("nDisplacedTracks", 2);
333  descriptions.add("l1tTrackJetsEmulator", desc);
334 }
335 
336 //define this as a plug-in
ap_uint< TTTrack_TrackWord::TrackBitWidths::kD0Size > d0_intern
std::vector< l1t::TkJetWord > TkJetWordCollection
Definition: TkJetWord.h:191
ap_int< kGlbEtaSize > glbeta_t
Definition: TkJetWord.h:59
const unsigned int kExtraGlobalPhiBit
const unsigned int DoubleToBit(double value, unsigned int maxBits, double step)
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
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
unsigned int eta_bin_firmwareStyle(int eta_word)
std::vector< TrackJetEmulationEtaPhiBin > clusters
ap_uint< kXtSize > nx_t
Definition: TkJetWord.h:63
unsigned int phi_bin_firmwareStyle(int phi_sector_raw, int phi_word)
const EDGetTokenT< L1TTTrackRefCollectionType > trackToken_
#define M_PI
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 b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(ConfigurationDescriptions &descriptions)
static constexpr double stepPhi0
HLT enums.
ap_uint< kNtSize > nt_t
Definition: TkJetWord.h:62
double a
Definition: hdecay.h:121
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
static constexpr double stepZ0
def move(src, dest)
Definition: eostools.py:511