CMS 3D CMS Logo

L1TrackJetEmulationProducer.cc
Go to the documentation of this file.
1 /* Software to emulate the hardware 2-layer jet-finding algorithm (fixed-point). *Layers 1 and 2*
2  *
3  * 2021
4  *
5  * Created based on Rishi Patel's L1TrackJetProducer.cc file.
6  * Authors: Samuel Edwin Leigh, Tyler Wu
7  * Rutgers, the State University of New Jersey
8  * Revolutionary for 250 years
9  */
10 
11 //Holds data from tracks, converted from their integer versions.
12 
13 // system include files
14 
23 // user include files
35 
37 
38 #include <memory>
39 #include <iostream>
40 #include <fstream>
41 #include <cmath>
42 #include <cstdlib>
43 #include <string>
44 #include <cstdlib>
45 #include "TH1D.h"
46 #include "TH2D.h"
47 #include <TMath.h>
48 #include <ap_int.h>
50 
51 using namespace std;
52 using namespace edm;
54 public:
55  explicit L1TrackJetEmulationProducer(const ParameterSet &);
56  ~L1TrackJetEmulationProducer() override;
58  typedef vector<L1TTTrackType> L1TTTrackCollectionType;
59 
60  static void fillDescriptions(ConfigurationDescriptions &descriptions);
61  bool trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2);
62  void L2_cluster(vector<Ptr<L1TTTrackType>> L1TrkPtrs_, TrackJetEmulationMaxZBin &mzb);
63  virtual TrackJetEmulationEtaPhiBin *L1_cluster(TrackJetEmulationEtaPhiBin *phislice);
64 
65 private:
66  void beginStream(StreamID) override;
67  void produce(Event &, const EventSetup &) override;
68  void endStream() override;
69 
70  // ----------member data ---------------------------
71 
72  vector<Ptr<L1TTTrackType>> L1TrkPtrs_;
73  vector<int> zBinCount_;
74  vector<int> tdtrk_;
75  const float MaxDzTrackPV;
76  const float trkZMax_;
77  const float trkPtMax_;
78  const float trkPtMin_;
79  const float trkEtaMax_;
80  const float trkChi2dofMax_;
81  const float trkBendChi2Max_;
82  const int trkNPSStubMin_;
83  const double minTrkJetpT_;
84  int etaBins_;
85  int phiBins_;
86  int zBins_;
93  bool displaced_;
99 
100  float PVz;
104 
108 };
109 
111  : MaxDzTrackPV((float)iConfig.getParameter<double>("MaxDzTrackPV")),
112  trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
113  trkPtMax_((float)iConfig.getParameter<double>("trk_ptMax")),
114  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
115  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
116  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
117  trkBendChi2Max_((float)iConfig.getParameter<double>("trk_bendChi2Max")),
118  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
119  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
120  etaBins_((int)iConfig.getParameter<int>("etaBins")),
121  phiBins_((int)iConfig.getParameter<int>("phiBins")),
122  zBins_((int)iConfig.getParameter<int>("zBins")),
123  d0CutNStubs4_((float)iConfig.getParameter<double>("d0_cutNStubs4")),
124  d0CutNStubs5_((float)iConfig.getParameter<double>("d0_cutNStubs5")),
125  lowpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
126  lowpTJetMinpT_((float)iConfig.getParameter<double>("lowpTJetMinpT")),
127  highpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
128  highpTJetMinpT_((float)iConfig.getParameter<double>("highpTJetMinpT")),
129  displaced_(iConfig.getParameter<bool>("displaced")),
130  nStubs4DisplacedChi2_((float)iConfig.getParameter<double>("nStubs4DisplacedChi2")),
131  nStubs5DisplacedChi2_((float)iConfig.getParameter<double>("nStubs5DisplacedChi2")),
132  nStubs4DisplacedBend_((float)iConfig.getParameter<double>("nStubs4Displacedbend")),
133  nStubs5DisplacedBend_((float)iConfig.getParameter<double>("nStubs5Displacedbend")),
134  nDisplacedTracks_((int)iConfig.getParameter<int>("nDisplacedTracks")),
135  trackToken_(consumes<vector<TTTrack<Ref_Phase2TrackerDigi_>>>(iConfig.getParameter<InputTag>("L1TrackInputTag"))),
136  PVtxToken_(consumes<l1t::VertexWordCollection>(iConfig.getParameter<InputTag>("VertexInputTag"))),
137  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
138  zStep_ = convert::makeZ0(2.0 * trkZMax_ / (zBins_ + 1));
139  etaStep_ = convert::makeGlbEta(2.0 * trkEtaMax_ / etaBins_); //etaStep is the width of an etabin
141 
142  if (displaced_)
143  produces<l1t::TkJetWordCollection>("L1TrackJetsExtended");
144  else
145  produces<l1t::TkJetWordCollection>("L1TrackJets");
146 }
147 
149 
151  unique_ptr<l1t::TkJetWordCollection> L1L1TrackJetProducer(new l1t::TkJetWordCollection);
152 
153  // For TTStubs
154  const TrackerTopology &tTopo = iSetup.getData(tTopoToken_);
155 
157  iEvent.getByToken(trackToken_, TTTrackHandle);
158  vector<TTTrack<Ref_Phase2TrackerDigi_>>::const_iterator iterL1Track;
159 
161  iEvent.getByToken(PVtxToken_, PVtx);
162  float PVz = (PVtx->at(0)).z0();
163 
164  L1TrkPtrs_.clear();
165  zBinCount_.clear();
166  tdtrk_.clear();
167  unsigned int this_l1track = 0;
168  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
169  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
170  this_l1track++;
171  float trk_pt = trkPtr->momentum().perp();
172  int trk_nstubs = (int)trkPtr->getStubRefs().size();
173  float trk_chi2dof = trkPtr->chi2Red();
174  float trk_d0 = trkPtr->d0();
175  float trk_bendchi2 = trkPtr->stubPtConsistency();
176  float trk_z0 = trkPtr->z0();
177 
178  int trk_nPS = 0;
179  for (int istub = 0; istub < trk_nstubs; istub++) { // loop over the stubs
180  DetId detId(trkPtr->getStubRefs().at(istub)->getDetId());
181  if (detId.det() == DetId::Detector::Tracker) {
182  if ((detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3) ||
183  (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9))
184  trk_nPS++;
185  }
186  }
187 
188  if (trk_nPS < trkNPSStubMin_)
189  continue;
190  if (!trackQualityCuts(trk_nstubs, trk_chi2dof, trk_bendchi2))
191  continue;
192  if (std::abs(trk_z0 - PVz) > MaxDzTrackPV)
193  continue;
194  if (std::abs(trk_z0) > trkZMax_)
195  continue;
196  if (std::abs(trkPtr->momentum().eta()) > trkEtaMax_)
197  continue;
198  if (trk_pt < trkPtMin_)
199  continue;
200  L1TrkPtrs_.push_back(trkPtr);
201  zBinCount_.push_back(0);
202 
203  if ((std::abs(trk_d0) > d0CutNStubs5_ && trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
204  (trk_nstubs == 4 && std::abs(trk_d0) > d0CutNStubs4_ && d0CutNStubs4_ >= 0))
205  tdtrk_.push_back(1); //displaced track
206  else
207  tdtrk_.push_back(0); // not displaced track
208  }
209 
210  if (!L1TrkPtrs_.empty()) {
212 
213  L2_cluster(L1TrkPtrs_, mzb);
214  if (mzb.clusters != nullptr) {
215  for (int j = 0; j < mzb.nclust; ++j) {
216  //FILL Two Layer Jets for Jet Collection
218  continue; //protects against reading bad memory
219  if (mzb.clusters[j].ntracks < 1)
220  continue;
221  if (mzb.clusters[j].ntracks > 5000)
222  continue;
227  l1t::TkJetWord::nt_t totalntracks_ = mzb.clusters[j].ntracks;
228  l1t::TkJetWord::nx_t totalxtracks_ = mzb.clusters[j].nxtracks;
229  l1t::TkJetWord::tkjetunassigned_t unassigned_ = 0;
230 
231  l1t::TkJetWord trkJet(jetPt, jetEta, jetPhi, jetZ0, totalntracks_, totalxtracks_, unassigned_);
232  //trkJet.setDispCounters(DispCounters);
233  L1L1TrackJetProducer->push_back(trkJet);
234  }
235  } else if (mzb.clusters == nullptr) {
236  edm::LogWarning("L1TrackJetEmulationProducer") << "mzb.clusters Not Assigned!\n";
237  }
238 
239  if (displaced_)
240  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
241  else
242  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
243  delete[] mzb.clusters;
244  } else if (L1TrkPtrs_.empty()) {
245  edm::LogWarning("L1TrackJetEmulationProducer") << "L1TrkPtrs Not Assigned!\n";
246  if (displaced_)
247  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
248  else
249  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
250  }
251 }
252 
254  enum TrackBitWidths {
255  kEtaSize = 16, // Width of z-position (40cm / 0.1)
256  kEtaMagSize = 3, // Width of z-position magnitude (signed)
257  kPtSize = 14, // Width of pt
258  kPtMagSize = 9, // Width of pt magnitude (unsigned)
259  kPhiSize = 12,
260  kPhiMagSize = 1,
261  };
262 
263  const int nz = zBins_ + 1;
264  TrackJetEmulationMaxZBin all_zBins[nz];
265  static TrackJetEmulationMaxZBin mzbtemp;
266  for (int z = 0; z < nz; ++z)
267  all_zBins[z] = mzbtemp;
268 
270  z0_intern zmax = zmin + 2 * zStep_;
271 
272  TrackJetEmulationEtaPhiBin epbins[phiBins_][etaBins_]; // create grid of phiBins
276  for (int i = 0; i < phiBins_; ++i) {
278  for (int j = 0; j < etaBins_; ++j) {
279  phimin = phi;
280  phimax = phi + phiStep_;
281  etamin = eta;
282  eta = eta + etaStep_;
283  etamax = eta;
284  epbins[i][j].phi = (phimin + phimax) / 2;
285  epbins[i][j].eta = (etamin + etamax) / 2;
286  epbins[i][j].pTtot = 0;
287  epbins[i][j].ntracks = 0;
288  epbins[i][j].nxtracks = 0;
289  } // for each etabin
290  phi = phi + phiStep_;
291  } // for each phibin (finished creating epbins)
292 
293  mzb = all_zBins[0];
294  mzb.ht = 0;
295  int ntracks = L1TrkPtrs_.size();
296  // uninitalized arrays
299 
300  for (int zbin = 0; zbin < zBins_; ++zbin) {
301  for (int i = 0; i < phiBins_; ++i) { //First initialize pT, numtracks, used to 0 (or false)
302  for (int j = 0; j < etaBins_; ++j) {
303  epbins[i][j].pTtot = 0;
304  epbins[i][j].used = false;
305  epbins[i][j].ntracks = 0;
306  epbins[i][j].nxtracks = 0;
307  } //for each etabin
308  L1clusters[i] = epbins[i];
309  } //for each phibin
310 
311  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
312  ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize, AP_RND_CONV, AP_SAT> inputTrkPt = 0;
313  inputTrkPt.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1,
314  TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
315  pt_intern trkpt = inputTrkPt;
316 
317  ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize, AP_RND_CONV, AP_SAT> trketainput = 0;
318  trketainput.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kTanlMSB,
319  TTTrack_TrackWord::TrackBitLocations::kTanlLSB);
320  ap_ufixed<64 + ETA_EXTRABITS, 32 + ETA_EXTRABITS> eta_conv =
321  1.0 / convert::ETA_LSB; //conversion factor from input eta format to internal format
322  glbeta_intern trketa = eta_conv * trketainput;
323 
324  int Sector = L1TrkPtrs_[k]->phiSector();
325  glbphi_intern trkphiSector = 0;
326  if (Sector < 5) {
327  trkphiSector = Sector * convert::makeGlbPhi(2.0 * M_PI / 9.0);
328  } else {
329  trkphiSector =
330  convert::makeGlbPhi(-1.0 * M_PI + M_PI / 9.0) + (Sector - 5) * convert::makeGlbPhi(2.0 * M_PI / 9.0);
331  }
332 
333  ap_fixed<TrackBitWidths::kPhiSize, TrackBitWidths::kPhiMagSize, AP_RND_CONV, AP_SAT> trkphiinput = 0;
334  trkphiinput.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kPhiMSB,
335  TTTrack_TrackWord::TrackBitLocations::kPhiLSB);
336  ap_ufixed<64 + PHI_EXTRABITS, 32 + PHI_EXTRABITS> phi_conv =
338  (1 << (TrackBitWidths::kPhiSize - TrackBitWidths::kPhiMagSize));
339  //phi_conv is a conversion factor from input phi format to the internal format
340  glbeta_intern localphi = phi_conv * trkphiinput;
341  glbeta_intern trkphi = localphi + trkphiSector;
342 
343  ap_int<TTTrack_TrackWord::TrackBitWidths::kZ0Size> inputTrkZ0 = L1TrkPtrs_[k]->getTrackWord()(
344  TTTrack_TrackWord::TrackBitLocations::kZ0MSB, TTTrack_TrackWord::TrackBitLocations::kZ0LSB);
345  ap_ufixed<32 + Z0_EXTRABITS, 1 + Z0_EXTRABITS> z0_conv =
346  TTTrack_TrackWord::stepZ0 / convert::Z0_LSB; //conversion factor from input z format to internal format
347  z0_intern trkZ = z0_conv * inputTrkZ0;
348 
349  for (int i = 0; i < phiBins_; ++i) {
350  for (int j = 0; j < etaBins_; ++j) {
351  L2cluster[k] = epbins[i][j];
352  if ((zmin <= trkZ && zmax >= trkZ) &&
353  ((epbins[i][j].eta - etaStep_ / 2 < trketa && epbins[i][j].eta + etaStep_ / 2 >= trketa) &&
354  epbins[i][j].phi - phiStep_ / 2 < trkphi && epbins[i][j].phi + phiStep_ / 2 >= trkphi &&
355  (zBinCount_[k] != 2))) {
356  zBinCount_.at(k) = zBinCount_.at(k) + 1;
357 
358  if (trkpt < convert::makePtFromFloat(trkPtMax_))
359  epbins[i][j].pTtot += trkpt;
360  else
362  ++epbins[i][j].ntracks;
363  //x-bit is currently not used in firmware, so we leave nxtracks = 0 for now
364  } // if right bin
365  } // for each phibin: j loop
366  } // for each phibin: i loop
367  } // end loop over tracks
368 
369  for (int phislice = 0; phislice < phiBins_; ++phislice) {
370  L1clusters[phislice] = L1_cluster(epbins[phislice]);
371  for (int ind = 0; L1clusters[phislice][ind].pTtot != 0; ++ind) {
372  L1clusters[phislice][ind].used = false;
373  }
374  }
375 
376  //Create clusters array to hold output cluster data for Layer2; can't have more clusters than tracks.
377  //Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used.
378  pt_intern hipT = 0;
379  int nclust = 0;
380  int phibin = 0;
381  int imax = -1;
382  int index1; //index of clusters array for each phislice
383  pt_intern E1 = 0;
384  pt_intern E0 = 0;
385  pt_intern E2 = 0;
386  l1t::TkJetWord::nt_t ntrk1, ntrk2;
387  l1t::TkJetWord::nx_t nxtrk1, nxtrk2;
388  int used1, used2, used3, used4;
389 
390  for (phibin = 0; phibin < phiBins_; ++phibin) { //Find eta-phibin with highest pT
391  while (true) {
392  hipT = 0;
393  for (index1 = 0; L1clusters[phibin][index1].pTtot > 0; ++index1) {
394  if (!L1clusters[phibin][index1].used && L1clusters[phibin][index1].pTtot >= hipT) {
395  hipT = L1clusters[phibin][index1].pTtot;
396  imax = index1;
397  }
398  } // for each index within the phibin
399 
400  if (hipT == 0)
401  break; // If highest pT is 0, all bins are used
402  E0 = hipT; // E0 is pT of first phibin of the cluster
403  E1 = 0;
404  E2 = 0;
405  ntrk1 = 0;
406  ntrk2 = 0;
407  nxtrk1 = 0;
408  nxtrk2 = 0;
409  L2cluster[nclust] = L1clusters[phibin][imax];
410 
411  L1clusters[phibin][imax].used = true;
412  // Add pT of upper neighbor
413  // E1 is pT of the middle phibin (should be highest pT)
414  if (phibin != phiBins_ - 1) {
415  used1 = -1;
416  used2 = -1;
417  for (index1 = 0; L1clusters[phibin + 1][index1].pTtot != 0; ++index1) {
418  if (L1clusters[phibin + 1][index1].used)
419  continue;
420  if (L1clusters[phibin + 1][index1].eta - L1clusters[phibin][imax].eta <= 3 * etaStep_ / 2 &&
421  L1clusters[phibin][imax].eta - L1clusters[phibin + 1][index1].eta <= 3 * etaStep_ / 2) {
422  E1 += L1clusters[phibin + 1][index1].pTtot;
423  ntrk1 += L1clusters[phibin + 1][index1].ntracks;
424  nxtrk1 += L1clusters[phibin + 1][index1].nxtracks;
425  if (used1 < 0)
426  used1 = index1;
427  else
428  used2 = index1;
429  } // if cluster is within one phibin
430  } // for each cluster in above phibin
431 
432  if (E1 < E0) { // if E1 isn't higher, E0 and E1 are their own cluster
433  L2cluster[nclust].pTtot += E1;
434  L2cluster[nclust].ntracks += ntrk1;
435  L2cluster[nclust].nxtracks += nxtrk1;
436  if (used1 >= 0)
437  L1clusters[phibin + 1][used1].used = true;
438  if (used2 >= 0)
439  L1clusters[phibin + 1][used2].used = true;
440  nclust++;
441  continue;
442  }
443 
444  if (phibin != phiBins_ - 2) { // E2 will be the pT of the third phibin (should be lower than E1)
445  used3 = -1;
446  used4 = -1;
447  for (index1 = 0; L1clusters[phibin + 2][index1].pTtot != 0; ++index1) {
448  if (L1clusters[phibin + 2][index1].used)
449  continue;
450  if (L1clusters[phibin + 2][index1].eta - L1clusters[phibin][imax].eta <= 3 * etaStep_ / 2 &&
451  L1clusters[phibin][imax].eta - L1clusters[phibin + 2][index1].eta <= 3 * etaStep_ / 2) {
452  E2 += L1clusters[phibin + 2][index1].pTtot;
453  ntrk2 += L1clusters[phibin + 2][index1].ntracks;
454  nxtrk2 += L1clusters[phibin + 2][index1].nxtracks;
455  if (used3 < 0)
456  used3 = index1;
457  else
458  used4 = index1;
459  }
460  }
461  // if indeed E2 < E1, add E1 and E2 to E0, they're all a cluster together
462  // otherwise, E0 is its own cluster
463  if (E2 < E1) {
464  L2cluster[nclust].pTtot += E1 + E2;
465  L2cluster[nclust].ntracks += ntrk1 + ntrk2;
466  L2cluster[nclust].nxtracks += nxtrk1 + nxtrk2;
467  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
468  if (used1 >= 0)
469  L1clusters[phibin + 1][used1].used = true;
470  if (used2 >= 0)
471  L1clusters[phibin + 1][used2].used = true;
472  if (used3 >= 0)
473  L1clusters[phibin + 2][used3].used = true;
474  if (used4 >= 0)
475  L1clusters[phibin + 2][used4].used = true;
476  }
477  nclust++;
478  continue;
479  } // end Not phiBins-2
480  else {
481  L2cluster[nclust].pTtot += E1;
482  L2cluster[nclust].ntracks += ntrk1;
483  L2cluster[nclust].nxtracks += nxtrk1;
484  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
485  if (used1 >= 0)
486  L1clusters[phibin + 1][used1].used = true;
487  if (used2 >= 0)
488  L1clusters[phibin + 1][used2].used = true;
489  nclust++;
490  continue;
491  }
492  } //End not last phibin(23)
493  else { //if it is phibin 23
494  L1clusters[phibin][imax].used = true;
495  nclust++;
496  }
497  } // while hipT not 0
498  } // for each phibin
499 
500  for (phibin = 0; phibin < phiBins_; ++phibin)
501  delete[] L1clusters[phibin];
502 
503  // Now merge clusters, if necessary
504  for (int m = 0; m < nclust - 1; ++m) {
505  for (int n = m + 1; n < nclust; ++n) {
506  if (L2cluster[n].eta == L2cluster[m].eta &&
507  ((L2cluster[n].phi - L2cluster[m].phi < 3 * phiStep_ / 2 &&
508  L2cluster[m].phi - L2cluster[n].phi < 3 * phiStep_ / 2) ||
509  (L2cluster[n].phi - L2cluster[m].phi > 2 * convert::makeGlbPhi(M_PI) - phiStep_ ||
510  L2cluster[m].phi - L2cluster[n].phi > 2 * convert::makeGlbPhi(M_PI) - phiStep_))) {
511  if (L2cluster[n].pTtot > L2cluster[m].pTtot) {
512  L2cluster[m].phi = L2cluster[n].phi;
513  }
514  L2cluster[m].pTtot += L2cluster[n].pTtot;
515  L2cluster[m].ntracks += L2cluster[n].ntracks;
516  L2cluster[m].nxtracks += L2cluster[n].nxtracks;
517  for (int m1 = n; m1 < nclust - 1; ++m1) {
518  L2cluster[m1] = L2cluster[m1 + 1];
519  }
520  nclust--;
521  m = -1;
522  break; //?????
523  } // end if clusters neighbor in eta
524  }
525  } // end for (m) loop
526  // sum up all pTs in this zbin to find ht
527 
528  pt_intern ht = 0;
529  for (int k = 0; k < nclust; ++k) {
530  if (L2cluster[k].pTtot > convert::makePtFromFloat(lowpTJetMinpT_) &&
532  continue;
533  if (L2cluster[k].pTtot > convert::makePtFromFloat(highpTJetMinpT_) &&
535  continue;
536  if (L2cluster[k].pTtot > convert::makePtFromFloat(minTrkJetpT_))
537  ht += L2cluster[k].pTtot;
538  }
539  // if ht is larger than previous max, this is the new vertex zbin
540  all_zBins[zbin].znum = zbin;
541  all_zBins[zbin].clusters = new TrackJetEmulationEtaPhiBin[nclust];
542  all_zBins[zbin].nclust = nclust;
543  all_zBins[zbin].zbincenter = (zmin + zmax) / 2;
544  for (int k = 0; k < nclust; ++k) {
545  all_zBins[zbin].clusters[k].phi = L2cluster[k].phi;
546  all_zBins[zbin].clusters[k].eta = L2cluster[k].eta;
547  all_zBins[zbin].clusters[k].pTtot = L2cluster[k].pTtot;
548  all_zBins[zbin].clusters[k].ntracks = L2cluster[k].ntracks;
549  all_zBins[zbin].clusters[k].nxtracks = L2cluster[k].nxtracks;
550  }
551  all_zBins[zbin].ht = ht;
552  if (ht >= mzb.ht) {
553  mzb = all_zBins[zbin];
554  mzb.zbincenter = (zmin + zmax) / 2;
555  }
556  // Prepare for next zbin!
557  zmin = zmin + zStep_;
558  zmax = zmax + zStep_;
559  } // for each zbin
560 
561  for (int zbin = 0; zbin < zBins_; ++zbin) {
562  if (zbin == mzb.znum) {
563  continue;
564  }
565  delete[] all_zBins[zbin].clusters;
566  }
567 }
568 
571  for (int etabin = 0; etabin < etaBins_ / 2; ++etabin) {
572  clusters[etabin].pTtot = 0;
573  clusters[etabin].ntracks = 0;
574  clusters[etabin].nxtracks = 0;
575  clusters[etabin].phi = 0;
576  clusters[etabin].eta = 0;
577  clusters[etabin].used = false;
578  }
579 
580  if (clusters == nullptr)
581  edm::LogWarning("L1TrackJetEmulationProducer") << "Clusters memory not assigned!\n";
582 
583  // Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used
584  pt_intern my_pt, left_pt, right_pt, right2pt;
585  int nclust = 0;
586  right2pt = 0;
587  for (int etabin = 0; etabin < etaBins_; ++etabin) {
588  // assign values for my pT and neighbors' pT
589  if (phislice[etabin].used)
590  continue;
591  my_pt = phislice[etabin].pTtot;
592  if (etabin > 0 && !phislice[etabin - 1].used) {
593  left_pt = phislice[etabin - 1].pTtot;
594  } else
595  left_pt = 0;
596  if (etabin < etaBins_ - 1 && !phislice[etabin + 1].used) {
597  right_pt = phislice[etabin + 1].pTtot;
598  if (etabin < etaBins_ - 2 && !phislice[etabin + 2].used) {
599  right2pt = phislice[etabin + 2].pTtot;
600  } else
601  right2pt = 0;
602  } else
603  right_pt = 0;
604 
605  // if I'm not a cluster, move on
606  if (my_pt < left_pt || my_pt <= right_pt) {
607  // if unused pT in the left neighbor, spit it out as a cluster
608  if (left_pt > 0) {
609  clusters[nclust] = phislice[etabin - 1];
610  phislice[etabin - 1].used = true;
611  nclust++;
612  }
613  continue;
614  }
615 
616  // I guess I'm a cluster-- should I use my right neighbor?
617  // Note: left neighbor will definitely be used because if it
618  // didn't belong to me it would have been used already
619  clusters[nclust] = phislice[etabin];
620  phislice[etabin].used = true;
621  if (left_pt > 0) {
622  clusters[nclust].pTtot += left_pt;
623  clusters[nclust].ntracks += phislice[etabin - 1].ntracks;
624  clusters[nclust].nxtracks += phislice[etabin - 1].nxtracks;
625  }
626  if (my_pt >= right2pt && right_pt > 0) {
627  clusters[nclust].pTtot += right_pt;
628  clusters[nclust].ntracks += phislice[etabin + 1].ntracks;
629  clusters[nclust].nxtracks += phislice[etabin + 1].nxtracks;
630  phislice[etabin + 1].used = true;
631  }
632  nclust++;
633  } // for each etabin
634 
635  // Now merge clusters, if necessary
636  for (int m = 0; m < nclust - 1; ++m) {
637  if (clusters[m + 1].eta - clusters[m].eta < 3 * etaStep_ / 2 &&
638  clusters[m].eta - clusters[m + 1].eta < 3 * etaStep_ / 2) {
639  if (clusters[m + 1].pTtot > clusters[m].pTtot) {
640  clusters[m].eta = clusters[m + 1].eta;
641  }
642  clusters[m].pTtot += clusters[m + 1].pTtot;
643  clusters[m].ntracks += clusters[m + 1].ntracks; // Previous version didn't add tracks when merging
644  clusters[m].nxtracks += clusters[m + 1].nxtracks;
645  for (int m1 = m + 1; m1 < nclust - 1; ++m1)
646  clusters[m1] = clusters[m1 + 1];
647  nclust--;
648  m = -1;
649  } // end if clusters neighbor in eta
650  } // end for (m) loop
651 
652  for (int i = nclust; i < etaBins_ / 2; ++i) // zero out remaining unused clusters
653  clusters[i].pTtot = 0;
654  return clusters;
655 }
656 
658 
660 
661 bool L1TrackJetEmulationProducer::trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2) {
662  bool PassQuality = false;
663  if (trk_bendchi2 < trkBendChi2Max_ && trk_chi2 < trkChi2dofMax_ && trk_nstub >= 4 && !displaced_)
664  PassQuality = true;
665  if (displaced_ && trk_bendchi2 < nStubs4DisplacedBend_ && trk_chi2 < nStubs4DisplacedChi2_ && trk_nstub == 4)
666  PassQuality = true;
667  if (displaced_ && trk_bendchi2 < nStubs5DisplacedBend_ && trk_chi2 < nStubs5DisplacedChi2_ && trk_nstub > 4)
668  PassQuality = true;
669  return PassQuality;
670 }
671 
673  //The following says we do not know what parameters are allowed so do no validation
674  // Please change this to state exactly what you do use, even if it is no parameters
676  desc.setUnknown();
677  descriptions.addDefault(desc);
678 }
679 
680 //define this as a plug-in
static const float ETA_LSB_POW
pt_intern makePtFromFloat(float pt)
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:156
static const float Z0_LSB_POW
T perp() const
Definition: PV3DBase.h:69
ap_int< kGlbEtaSize > glbeta_t
Definition: TkJetWord.h:55
ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
glbeta_intern makeGlbEta(float eta)
ap_int< 14+PHI_EXTRABITS > glbphi_intern
T eta() const
Definition: PV3DBase.h:73
static const float ETA_LSB
z0_intern makeZ0(float z0)
vector< L1TTTrackType > L1TTTrackCollectionType
ap_int< kGlbPhiSize > glbphi_t
Definition: TkJetWord.h:56
delete x;
Definition: CaloConfig.h:22
void L2_cluster(vector< Ptr< L1TTTrackType >> L1TrkPtrs_, TrackJetEmulationMaxZBin &mzb)
bool trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2)
ap_ufixed< kPtSize, kPtMagSize, AP_TRN, AP_SAT > pt_t
Definition: TkJetWord.h:54
static const float PHI_LSB_POW
ap_uint< TkJetBitWidths::kUnassignedSize > tkjetunassigned_t
Definition: TkJetWord.h:60
ap_int< kZ0Size > z0_t
Definition: TkJetWord.h:57
static void fillDescriptions(ConfigurationDescriptions &descriptions)
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
ap_ufixed< 16+PT_EXTRABITS, 11, AP_TRN, AP_SAT > pt_intern
ap_int< 10+Z0_EXTRABITS > z0_intern
ap_int< 14+ETA_EXTRABITS > glbeta_intern
glbphi_intern makeGlbPhi(float phi)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
static const float PHI_LSB
void produce(Event &, const EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< VertexWord > VertexWordCollection
Definition: VertexWord.h:195
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
const EDGetTokenT< l1t::VertexWordCollection > PVtxToken_
ap_uint< kXtSize > nx_t
Definition: TkJetWord.h:59
std::vector< edm::Ref< edmNew::DetSetVector< TTStub< T > >, TTStub< T > > > getStubRefs() const
Track components.
Definition: TTTrack.h:93
#define M_PI
Definition: DetId.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
double stubPtConsistency() const
StubPtConsistency.
Definition: TTTrack.h:417
static constexpr double stepPhi0
double d0() const
Track d0.
Definition: TTTrack.h:325
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
HLT enums.
double chi2Red() const
Chi2 reduced.
Definition: TTTrack.h:359
ap_uint< kNtSize > nt_t
Definition: TkJetWord.h:58
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
double z0() const
Track z0.
Definition: TTTrack.h:330
TrackJetEmulationEtaPhiBin * clusters
unsigned int tidRing(const DetId &id) const
const EDGetTokenT< vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Log< level::Warning, false > LogWarning
static constexpr double stepZ0
static constexpr auto TID
static const float Z0_LSB
def move(src, dest)
Definition: eostools.py:511
L1TrackJetEmulationProducer(const ParameterSet &)
virtual TrackJetEmulationEtaPhiBin * L1_cluster(TrackJetEmulationEtaPhiBin *phislice)