CMS 3D CMS Logo

L1TrackJetProducer.cc
Go to the documentation of this file.
1 // Original Author: Rishi Patel
2 // Modifications: George Karathanasis, georgios.karathanasis@cern.ch, CU Boulder
3 // Created: Wed, 01 Aug 2018 14:01:41 GMT
4 //
5 // Track jets are clustered in a two-layer process, first by clustering in phi,
6 // then by clustering in eta
7 // Introduction to object (p10-13):
8 // https://indico.cern.ch/event/791517/contributions/3341650/attachments/1818736/2973771/TrackBasedAlgos_L1TMadrid_MacDonald.pdf
9 
10 // system include files
11 
17 // user include files
29 
32 #include "TH1D.h"
33 #include "TH2D.h"
34 #include <TMath.h>
35 #include <cmath>
36 #include <cstdlib>
37 #include <fstream>
38 #include <iostream>
39 #include <memory>
40 #include <string>
41 
42 using namespace std;
43 using namespace edm;
44 using namespace l1t;
46 public:
47  explicit L1TrackJetProducer(const ParameterSet &);
48  ~L1TrackJetProducer() override;
50  typedef vector<L1TTTrackType> L1TTTrackCollectionType;
51 
52  static void fillDescriptions(ConfigurationDescriptions &descriptions);
53  bool trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2);
54  void L2_cluster(vector<Ptr<L1TTTrackType> > L1TrkPtrs_, vector<int> tdtrk_, MaxZBin &mzb);
55  virtual EtaPhiBin *L1_cluster(EtaPhiBin *phislice);
56 
57 private:
58  void beginStream(StreamID) override;
59  void produce(Event &, const EventSetup &) override;
60  void endStream() override;
61 
62  // ----------member data ---------------------------
63 
64  vector<Ptr<L1TTTrackType> > L1TrkPtrs_;
65  vector<int> zBinCount_;
66  vector<int> tdtrk_;
67  const float trkZMax_;
68  const float trkPtMax_;
69  const float trkPtMin_;
70  const float trkEtaMax_;
71  const float trkChi2dofMax_;
72  const float trkBendChi2Max_;
73  const int trkNPSStubMin_;
74  const double minTrkJetpT_;
75  const double minJetEtLowPt_;
76  const double minJetEtHighPt_;
77  int etaBins_;
78  int phiBins_;
79  int zBins_;
84  bool displaced_;
90  float dzPVTrk_;
91 
95 
96  float zStep_;
97  float etaStep_;
98  float phiStep_;
99 };
100 
102  : trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
103  trkPtMax_((float)iConfig.getParameter<double>("trk_ptMax")),
104  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
105  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
106  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
107  trkBendChi2Max_((float)iConfig.getParameter<double>("trk_bendChi2Max")),
108  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
109  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
110  minJetEtLowPt_(iConfig.getParameter<double>("minJetEtLowPt")),
111  minJetEtHighPt_(iConfig.getParameter<double>("minJetEtHighPt")),
112  etaBins_((int)iConfig.getParameter<int>("etaBins")),
113  phiBins_((int)iConfig.getParameter<int>("phiBins")),
114  zBins_((int)iConfig.getParameter<int>("zBins")),
115  d0CutNStubs4_((float)iConfig.getParameter<double>("d0_cutNStubs4")),
116  d0CutNStubs5_((float)iConfig.getParameter<double>("d0_cutNStubs5")),
117  lowpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
118  highpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
119  displaced_(iConfig.getParameter<bool>("displaced")),
120  nStubs4DisplacedChi2_((float)iConfig.getParameter<double>("nStubs4DisplacedChi2")),
121  nStubs5DisplacedChi2_((float)iConfig.getParameter<double>("nStubs5DisplacedChi2")),
122  nStubs4DisplacedBend_((float)iConfig.getParameter<double>("nStubs4Displacedbend")),
123  nStubs5DisplacedBend_((float)iConfig.getParameter<double>("nStubs5Displacedbend")),
124  nDisplacedTracks_((int)iConfig.getParameter<int>("nDisplacedTracks")),
125  dzPVTrk_((float)iConfig.getParameter<double>("MaxDzTrackPV")),
126  trackToken_(
127  consumes<vector<TTTrack<Ref_Phase2TrackerDigi_> > >(iConfig.getParameter<InputTag>("L1TrackInputTag"))),
128  PVtxToken_(consumes<vector<l1t::Vertex> >(iConfig.getParameter<InputTag>("L1PVertexCollection"))),
129  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
130  zStep_ = 2.0 * trkZMax_ / (zBins_ + 1);
131  etaStep_ = 2.0 * trkEtaMax_ / etaBins_; //etaStep is the width of an etabin
132  phiStep_ = 2 * M_PI / phiBins_;
133 
134  if (displaced_)
135  produces<TkJetCollection>("L1TrackJetsExtended");
136  else
137  produces<TkJetCollection>("L1TrackJets");
138 }
139 
141 
143  unique_ptr<TkJetCollection> L1L1TrackJetProducer(new TkJetCollection);
144 
145  // For TTStubs
146  const TrackerTopology &tTopo = iSetup.getData(tTopoToken_);
147 
149  iEvent.getByToken(trackToken_, TTTrackHandle);
150  vector<TTTrack<Ref_Phase2TrackerDigi_> >::const_iterator iterL1Track;
151 
153  iEvent.getByToken(PVtxToken_, PVtx);
154  float PVz = (PVtx->at(0)).z0();
155 
156  L1TrkPtrs_.clear();
157  zBinCount_.clear();
158  tdtrk_.clear();
159  unsigned int this_l1track = 0;
160  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
161  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
162  this_l1track++;
163  float trk_pt = trkPtr->momentum().perp();
164  int trk_nstubs = (int)trkPtr->getStubRefs().size();
165  float trk_chi2dof = trkPtr->chi2Red();
166  float trk_d0 = trkPtr->d0();
167  float trk_bendchi2 = trkPtr->stubPtConsistency();
168  float trk_z0 = trkPtr->z0();
169 
170  int trk_nPS = 0;
171  for (int istub = 0; istub < trk_nstubs; istub++) { // loop over the stubs
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  if (trk_nPS < trkNPSStubMin_)
180  continue;
181  if (!trackQualityCuts(trk_nstubs, trk_chi2dof, trk_bendchi2))
182  continue;
183  if (std::abs(PVz - trk_z0) > dzPVTrk_)
184  continue;
185  if (std::abs(trk_z0) > trkZMax_)
186  continue;
187  if (std::abs(trkPtr->momentum().eta()) > trkEtaMax_)
188  continue;
189  if (trk_pt < trkPtMin_)
190  continue;
191  L1TrkPtrs_.push_back(trkPtr);
192  zBinCount_.push_back(0);
193 
194  if ((std::abs(trk_d0) > d0CutNStubs5_ && trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
195  (trk_nstubs == 4 && std::abs(trk_d0) > d0CutNStubs4_ && d0CutNStubs4_ >= 0))
196  tdtrk_.push_back(1); //displaced track
197  else
198  tdtrk_.push_back(0); // not displaced track
199  }
200 
201  if (!L1TrkPtrs_.empty()) {
202  MaxZBin mzb;
203 
205  vector<Ptr<L1TTTrackType> > L1TrackAssocJet;
206  if (mzb.clusters != nullptr) {
207  for (int j = 0; j < mzb.nclust; ++j) {
208  //FILL Two Layer Jets for Jet Collection
209  if (mzb.clusters[j].pTtot <= trkPtMin_)
210  continue; //protects against reading bad memory
211  if (mzb.clusters[j].numtracks < 1)
212  continue;
213  if (mzb.clusters[j].numtracks > 5000)
214  continue;
215  float jetEta = mzb.clusters[j].eta;
216  float jetPhi = mzb.clusters[j].phi;
217  float jetPt = mzb.clusters[j].pTtot;
218  float jetPx = jetPt * cos(jetPhi);
219  float jetPy = jetPt * sin(jetPhi);
220  float jetPz = jetPt * sinh(jetEta);
221  float jetP = jetPt * cosh(jetEta);
222  int totalDisptrk = mzb.clusters[j].numtdtrks;
223  bool isDispJet = totalDisptrk >= nDisplacedTracks_;
224  math::XYZTLorentzVector jetP4(jetPx, jetPy, jetPz, jetP);
225  L1TrackAssocJet.clear();
226  for (unsigned int itrk = 0; itrk < mzb.clusters[j].trackidx.size(); itrk++) {
227  L1TrackAssocJet.push_back(L1TrkPtrs_[mzb.clusters[j].trackidx[itrk]]);
228  }
229 
230  TkJet trkJet(jetP4, L1TrackAssocJet, mzb.zbincenter, mzb.clusters[j].numtracks, 0, totalDisptrk, 0, isDispJet);
231  if (!L1TrackAssocJet.empty()) {
232  L1L1TrackJetProducer->push_back(trkJet);
233  }
234  }
235  }
236 
237  if (displaced_)
238  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
239  else
240  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
241  delete[] mzb.clusters;
242  }
243 }
244 
245 void L1TrackJetProducer::L2_cluster(vector<Ptr<L1TTTrackType> > L1TrkPtrs_, vector<int> tdtrk_, MaxZBin &mzb) {
246  const int nz = zBins_ + 1;
247  MaxZBin all_zBins[nz];
248  MaxZBin mzbtemp = {};
249  for (int z = 0; z < nz; ++z)
250  all_zBins[z] = mzbtemp;
251 
252  float zmin = -1.0 * trkZMax_;
253  float zmax = zmin + 2 * zStep_;
254 
255  EtaPhiBin epbins[phiBins_][etaBins_]; // create grid of phiBins
256  float phi = -1.0 * M_PI;
257  float eta;
258  float etamin, etamax, phimin, phimax;
259 
260  for (int i = 0; i < phiBins_; ++i) {
261  eta = -1.0 * trkEtaMax_;
262  for (int j = 0; j < etaBins_; ++j) {
263  phimin = phi;
264  phimax = phi + phiStep_;
265  etamin = eta;
266  eta = eta + etaStep_;
267  etamax = eta;
268  epbins[i][j].phi = (phimin + phimax) / 2.0;
269  epbins[i][j].eta = (etamin + etamax) / 2.0;
270  } // for each etabin
271  phi = phi + phiStep_;
272  } // for each phibin (finished creating epbins)
273 
274  mzb = all_zBins[0];
275  int ntracks = L1TrkPtrs_.size();
276 
277  // uninitalized arrays
278  EtaPhiBin *L1clusters[phiBins_];
279  EtaPhiBin L2cluster[ntracks];
280 
281  for (int zbin = 0; zbin < zBins_; ++zbin) {
282  for (int i = 0; i < phiBins_; ++i) { //First initialize pT, numtracks, used to 0 (or false)
283  for (int j = 0; j < etaBins_; ++j) {
284  epbins[i][j].pTtot = 0;
285  epbins[i][j].used = false;
286  epbins[i][j].numtracks = 0;
287  epbins[i][j].numttrks = 0;
288  epbins[i][j].numtdtrks = 0;
289  epbins[i][j].numttdtrks = 0;
290  epbins[i][j].trackidx.clear();
291  } //for each etabin
292  L1clusters[i] = epbins[i];
293  } //for each phibin
294 
295  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
296  float trkpt = L1TrkPtrs_[k]->momentum().perp();
297  float trketa = L1TrkPtrs_[k]->momentum().eta();
298  float trkphi = L1TrkPtrs_[k]->momentum().phi();
299  float trkZ = L1TrkPtrs_[k]->z0();
300 
301  for (int i = 0; i < phiBins_; ++i) {
302  for (int j = 0; j < etaBins_; ++j) {
303  L2cluster[k] = epbins[i][j];
304  if ((zmin <= trkZ && zmax >= trkZ) &&
305  ((epbins[i][j].eta - etaStep_ / 2.0 < trketa && epbins[i][j].eta + etaStep_ / 2.0 >= trketa) &&
306  epbins[i][j].phi - phiStep_ / 2.0 < trkphi && epbins[i][j].phi + phiStep_ / 2.0 >= trkphi &&
307  (zBinCount_[k] != 2))) {
308  zBinCount_.at(k) = zBinCount_.at(k) + 1;
309  if (trkpt < trkPtMax_)
310  epbins[i][j].pTtot += trkpt;
311  else
312  epbins[i][j].pTtot += trkPtMax_;
313  epbins[i][j].numtdtrks += tdtrk_[k];
314  epbins[i][j].trackidx.push_back(k);
315  ++epbins[i][j].numtracks;
316  } // if right bin
317  } // for each phibin: j loop
318  } // for each phibin: i loop
319  } // end loop over tracks
320 
321  for (int phislice = 0; phislice < phiBins_; ++phislice) {
322  L1clusters[phislice] = L1_cluster(epbins[phislice]);
323  if (L1clusters[phislice] != nullptr) {
324  for (int ind = 0; L1clusters[phislice][ind].pTtot != 0; ++ind) {
325  L1clusters[phislice][ind].used = false;
326  }
327  }
328  }
329 
330  //Create clusters array to hold output cluster data for Layer2; can't have more clusters than tracks.
331  //Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used.
332  float hipT = 0;
333  int nclust = 0;
334  int phibin = 0;
335  int imax = -1;
336  int index1; //index of clusters array for each phislice
337  float E1 = 0;
338  float E0 = 0;
339  float E2 = 0;
340  int trx1, trx2;
341  int tdtrk1, tdtrk2;
342  int used1, used2, used3, used4;
343 
344  for (phibin = 0; phibin < phiBins_; ++phibin) { //Find eta-phibin with highest pT
345  while (true) {
346  hipT = 0;
347  for (index1 = 0; L1clusters[phibin][index1].pTtot > 0; ++index1) {
348  if (!L1clusters[phibin][index1].used && L1clusters[phibin][index1].pTtot >= hipT) {
349  hipT = L1clusters[phibin][index1].pTtot;
350  imax = index1;
351  }
352  } // for each index within the phibin
353 
354  if (hipT == 0)
355  break; // If highest pT is 0, all bins are used
356  E0 = hipT; // E0 is pT of first phibin of the cluster
357  E1 = 0;
358  E2 = 0;
359  trx1 = 0;
360  trx2 = 0;
361  tdtrk1 = 0;
362  tdtrk2 = 0;
363  std::vector<unsigned int> trkidx1;
364  std::vector<unsigned int> trkidx2;
365  L2cluster[nclust] = L1clusters[phibin][imax];
366  L1clusters[phibin][imax].used = true;
367  // Add pT of upper neighbor
368  // E1 is pT of the middle phibin (should be highest pT)
369  if (phibin != phiBins_ - 1) {
370  used1 = -1;
371  used2 = -1;
372  for (index1 = 0; L1clusters[phibin + 1][index1].pTtot != 0; ++index1) {
373  if (L1clusters[phibin + 1][index1].used)
374  continue;
375 
376  if (std::abs(L1clusters[phibin + 1][index1].eta - L1clusters[phibin][imax].eta) <= 1.5 * etaStep_) {
377  E1 += L1clusters[phibin + 1][index1].pTtot;
378  trx1 += L1clusters[phibin + 1][index1].numtracks;
379  tdtrk1 += L1clusters[phibin + 1][index1].numtdtrks;
380  for (unsigned int itrk = 0; itrk < L1clusters[phibin + 1][index1].trackidx.size(); itrk++) {
381  trkidx1.push_back(L1clusters[phibin + 1][index1].trackidx[itrk]);
382  }
383  if (used1 < 0)
384  used1 = index1;
385  else
386  used2 = index1;
387  } // if cluster is within one phibin
388  } // for each cluster in above phibin
389 
390  if (E1 < E0) { // if E1 isn't higher, E0 and E1 are their own cluster
391  L2cluster[nclust].pTtot += E1;
392  L2cluster[nclust].numtracks += trx1;
393  L2cluster[nclust].numtdtrks += tdtrk1;
394  for (unsigned int itrk = 0; itrk < trkidx1.size(); itrk++) {
395  L2cluster[nclust].trackidx.push_back(trkidx1[itrk]);
396  }
397 
398  if (used1 >= 0)
399  L1clusters[phibin + 1][used1].used = true;
400  if (used2 >= 0)
401  L1clusters[phibin + 1][used2].used = true;
402  nclust++;
403  continue;
404  }
405 
406  if (phibin != phiBins_ - 2) { // E2 will be the pT of the third phibin (should be lower than E1)
407  used3 = -1;
408  used4 = -1;
409  for (index1 = 0; L1clusters[phibin + 2][index1].pTtot != 0; ++index1) {
410  if (L1clusters[phibin + 2][index1].used)
411  continue;
412  if (std::abs(L1clusters[phibin + 2][index1].eta - L1clusters[phibin][imax].eta) <= 1.5 * etaStep_) {
413  E2 += L1clusters[phibin + 2][index1].pTtot;
414  trx2 += L1clusters[phibin + 2][index1].numtracks;
415  tdtrk2 += L1clusters[phibin + 2][index1].numtdtrks;
416  for (unsigned int itrk = 0; itrk < L1clusters[phibin + 2][index1].trackidx.size(); itrk++)
417  trkidx2.push_back(L1clusters[phibin + 2][index1].trackidx[itrk]);
418 
419  if (used3 < 0)
420  used3 = index1;
421  else
422  used4 = index1;
423  }
424  }
425 
426  // if indeed E2 < E1, add E1 and E2 to E0, they're all a cluster together
427  // otherwise, E0 is its own cluster
428  if (E2 < E1) {
429  L2cluster[nclust].pTtot += E1 + E2;
430  L2cluster[nclust].numtracks += trx1 + trx2;
431  L2cluster[nclust].numtdtrks += tdtrk1 + tdtrk2;
432  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
433  for (unsigned int itrk = 0; itrk < trkidx1.size(); itrk++) {
434  L2cluster[nclust].trackidx.push_back(trkidx1[itrk]);
435  }
436  for (unsigned int itrk = 0; itrk < trkidx2.size(); itrk++) {
437  L2cluster[nclust].trackidx.push_back(trkidx2[itrk]);
438  }
439 
440  if (used1 >= 0)
441  L1clusters[phibin + 1][used1].used = true;
442  if (used2 >= 0)
443  L1clusters[phibin + 1][used2].used = true;
444  if (used3 >= 0)
445  L1clusters[phibin + 2][used3].used = true;
446  if (used4 >= 0)
447  L1clusters[phibin + 2][used4].used = true;
448  }
449  nclust++;
450  continue;
451  } // end Not phiBins-2
452  else {
453  L2cluster[nclust].pTtot += E1;
454  L2cluster[nclust].numtracks += trx1;
455  L2cluster[nclust].numtdtrks += tdtrk1;
456  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
457  for (unsigned int itrk = 0; itrk < trkidx1.size(); itrk++) {
458  L2cluster[nclust].trackidx.push_back(trkidx1[itrk]);
459  }
460 
461  if (used1 >= 0)
462  L1clusters[phibin + 1][used1].used = true;
463  if (used2 >= 0)
464  L1clusters[phibin + 1][used2].used = true;
465  nclust++;
466  continue;
467  }
468  } //End not last phibin(23)
469  else { //if it is phibin 23
470  L1clusters[phibin][imax].used = true;
471  nclust++;
472  }
473  } // while hipT not 0
474  } // for each phibin
475 
476  for (phibin = 0; phibin < phiBins_; ++phibin)
477  delete[] L1clusters[phibin];
478 
479  // Now merge clusters, if necessary
480  for (int m = 0; m < nclust - 1; ++m) {
481  for (int n = m + 1; n < nclust; ++n) {
482  if (L2cluster[n].eta == L2cluster[m].eta && (std::abs(L2cluster[n].phi - L2cluster[m].phi) < 1.5 * phiStep_ ||
483  std::abs(L2cluster[n].phi - L2cluster[m].phi) > 6.0)) {
484  if (L2cluster[n].pTtot > L2cluster[m].pTtot) {
485  L2cluster[m].phi = L2cluster[n].phi;
486  }
487  L2cluster[m].pTtot += L2cluster[n].pTtot;
488  L2cluster[m].numtracks += L2cluster[n].numtracks;
489  L2cluster[m].numtdtrks += L2cluster[n].numtdtrks;
490  for (unsigned int itrk = 0; itrk < L2cluster[n].trackidx.size(); itrk++)
491  L2cluster[m].trackidx.push_back(L2cluster[n].trackidx[itrk]);
492 
493  for (int m1 = n; m1 < nclust - 1; ++m1) {
494  L2cluster[m1] = L2cluster[m1 + 1];
495  }
496  nclust--;
497  m = -1;
498  break; //?????
499  } // end if clusters neighbor in eta
500  }
501  } // end for (m) loop
502 
503  // sum up all pTs in this zbin to find ht
504  // Require jet to have at least 2 tracks with a summed eT>50 GeV, or 3 tracks with summed eT>100 GeV
505  // in order to count toward HT
506  float ht = 0;
507  for (int k = 0; k < nclust; ++k) {
508  if (L2cluster[k].pTtot > minJetEtLowPt_ && L2cluster[k].numtracks < lowpTJetMinTrackMultiplicity_)
509  continue;
510  if (L2cluster[k].pTtot > minJetEtHighPt_ && L2cluster[k].numtracks < highpTJetMinTrackMultiplicity_)
511  continue;
512  if (L2cluster[k].pTtot > minTrkJetpT_)
513  ht += L2cluster[k].pTtot;
514  }
515 
516  // if ht is larger than previous max, this is the new vertex zbin
517  all_zBins[zbin].znum = zbin;
518  all_zBins[zbin].clusters = new EtaPhiBin[nclust];
519  all_zBins[zbin].nclust = nclust;
520  all_zBins[zbin].zbincenter = (zmin + zmax) / 2.0;
521  for (int k = 0; k < nclust; ++k) {
522  all_zBins[zbin].clusters[k].phi = L2cluster[k].phi;
523  all_zBins[zbin].clusters[k].eta = L2cluster[k].eta;
524  all_zBins[zbin].clusters[k].pTtot = L2cluster[k].pTtot;
525  all_zBins[zbin].clusters[k].numtracks = L2cluster[k].numtracks;
526  all_zBins[zbin].clusters[k].numtdtrks = L2cluster[k].numtdtrks;
527  for (unsigned int itrk = 0; itrk < L2cluster[k].trackidx.size(); itrk++)
528  all_zBins[zbin].clusters[k].trackidx.push_back(L2cluster[k].trackidx[itrk]);
529  }
530 
531  all_zBins[zbin].ht = ht;
532  if (ht >= mzb.ht) {
533  mzb = all_zBins[zbin];
534  mzb.zbincenter = (zmin + zmax) / 2.0;
535  }
536  // Prepare for next zbin!
537  zmin = zmin + zStep_;
538  zmax = zmax + zStep_;
539  } // for each zbin
540  for (int zbin = 0; zbin < zBins_; ++zbin) {
541  if (zbin == mzb.znum)
542  continue;
543  delete[] all_zBins[zbin].clusters;
544  }
545 }
546 
548  EtaPhiBin *clusters = new EtaPhiBin[etaBins_ / 2];
549  if (clusters == nullptr)
550  edm::LogWarning("L1TrackJetProducer") << "Clusters memory not assigned!\n";
551 
552  // Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used
553  float my_pt, left_pt, right_pt, right2pt;
554  int nclust = 0;
555  right2pt = 0;
556  for (int etabin = 0; etabin < etaBins_; ++etabin) {
557  // assign values for my pT and neighbors' pT
558  if (phislice[etabin].used)
559  continue;
560  my_pt = phislice[etabin].pTtot;
561  if (etabin > 0 && !phislice[etabin - 1].used) {
562  left_pt = phislice[etabin - 1].pTtot;
563  } else
564  left_pt = 0;
565  if (etabin < etaBins_ - 1 && !phislice[etabin + 1].used) {
566  right_pt = phislice[etabin + 1].pTtot;
567  if (etabin < etaBins_ - 2 && !phislice[etabin + 2].used) {
568  right2pt = phislice[etabin + 2].pTtot;
569  } else
570  right2pt = 0;
571  } else
572  right_pt = 0;
573 
574  // if I'm not a cluster, move on
575  if (my_pt < left_pt || my_pt <= right_pt) {
576  // if unused pT in the left neighbor, spit it out as a cluster
577  if (left_pt > 0) {
578  clusters[nclust] = phislice[etabin - 1];
579  phislice[etabin - 1].used = true;
580  nclust++;
581  }
582  continue;
583  }
584 
585  // I guess I'm a cluster-- should I use my right neighbor?
586  // Note: left neighbor will definitely be used because if it
587  // didn't belong to me it would have been used already
588  clusters[nclust] = phislice[etabin];
589  phislice[etabin].used = true;
590  if (left_pt > 0) {
591  if (clusters != nullptr) {
592  clusters[nclust].pTtot += left_pt;
593  clusters[nclust].numtracks += phislice[etabin - 1].numtracks;
594  clusters[nclust].numtdtrks += phislice[etabin - 1].numtdtrks;
595  for (unsigned int itrk = 0; itrk < phislice[etabin - 1].trackidx.size(); itrk++)
596  clusters[nclust].trackidx.push_back(phislice[etabin - 1].trackidx[itrk]);
597  }
598  }
599  if (my_pt >= right2pt && right_pt > 0) {
600  if (clusters != nullptr) {
601  clusters[nclust].pTtot += right_pt;
602  clusters[nclust].numtracks += phislice[etabin + 1].numtracks;
603  clusters[nclust].numtdtrks += phislice[etabin + 1].numtdtrks;
604  for (unsigned int itrk = 0; itrk < phislice[etabin + 1].trackidx.size(); itrk++)
605  clusters[nclust].trackidx.push_back(phislice[etabin + 1].trackidx[itrk]);
606  phislice[etabin + 1].used = true;
607  }
608  }
609  nclust++;
610  } // for each etabin
611 
612  // Now merge clusters, if necessary
613  for (int m = 0; m < nclust - 1; ++m) {
614  if (std::abs(clusters[m + 1].eta - clusters[m].eta) < 1.5 * etaStep_) {
615  if (clusters[m + 1].pTtot > clusters[m].pTtot) {
616  clusters[m].eta = clusters[m + 1].eta;
617  }
618  clusters[m].pTtot += clusters[m + 1].pTtot;
619  clusters[m].numtracks += clusters[m + 1].numtracks; // Previous version didn't add tracks when merging
620  clusters[m].numtdtrks += clusters[m + 1].numtdtrks;
621  for (unsigned int itrk = 0; itrk < clusters[m + 1].trackidx.size(); itrk++)
622  clusters[m].trackidx.push_back(clusters[m + 1].trackidx[itrk]);
623 
624  for (int m1 = m + 1; m1 < nclust - 1; ++m1)
625  clusters[m1] = clusters[m1 + 1];
626  nclust--;
627  m = -1;
628  } // end if clusters neighbor in eta
629  } // end for (m) loop
630 
631  for (int i = nclust; i < etaBins_ / 2; ++i) // zero out remaining unused clusters
632  clusters[i].pTtot = 0;
633  return clusters;
634 }
635 
637 
639 
640 bool L1TrackJetProducer::trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2) {
641  bool PassQuality = false;
642  if (trk_bendchi2 < trkBendChi2Max_ && trk_chi2 < trkChi2dofMax_ && trk_nstub >= 4 && !displaced_)
643  PassQuality = true;
644  if (displaced_ && trk_bendchi2 < nStubs4DisplacedBend_ && trk_chi2 < nStubs4DisplacedChi2_ && trk_nstub == 4)
645  PassQuality = true;
646  if (displaced_ && trk_bendchi2 < nStubs5DisplacedBend_ && trk_chi2 < nStubs5DisplacedChi2_ && trk_nstub > 4)
647  PassQuality = true;
648  return PassQuality;
649 }
650 
652  //The following says we do not know what parameters are allowed so do no validation
653  // Please change this to state exactly what you do use, even if it is no parameters
655  desc.setUnknown();
656  descriptions.addDefault(desc);
657 }
658 
659 //define this as a plug-in
bool trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2)
const int highpTJetMinTrackMultiplicity_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int tobLayer(const DetId &id) const
T perp() const
Definition: PV3DBase.h:69
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
Definition: PV3DBase.h:73
delete x;
Definition: CaloConfig.h:22
L1TrackJetProducer(const ParameterSet &)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
std::vector< unsigned int > trackidx
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
EtaPhiBin * clusters
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
vector< L1TTTrackType > L1TTTrackCollectionType
void addDefault(ParameterSetDescription const &psetDescription)
static void fillDescriptions(ConfigurationDescriptions &descriptions)
const double minJetEtLowPt_
float zbincenter
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual EtaPhiBin * L1_cluster(EtaPhiBin *phislice)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
vector< int > zBinCount_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void endStream() override
static constexpr auto TOB
std::vector< edm::Ref< edmNew::DetSetVector< TTStub< T > >, TTStub< T > > > getStubRefs() const
Track components.
Definition: TTTrack.h:93
#define M_PI
void beginStream(StreamID) override
const edm::EDGetTokenT< std::vector< l1t::Vertex > > PVtxToken_
Definition: DetId.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(Event &, const EventSetup &) override
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
double stubPtConsistency() const
StubPtConsistency.
Definition: TTTrack.h:417
double d0() const
Track d0.
Definition: TTTrack.h:325
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
const double minJetEtHighPt_
const int lowpTJetMinTrackMultiplicity_
const EDGetTokenT< vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
HLT enums.
double chi2Red() const
Chi2 reduced.
Definition: TTTrack.h:359
void L2_cluster(vector< Ptr< L1TTTrackType > > L1TrkPtrs_, vector< int > tdtrk_, MaxZBin &mzb)
double z0() const
Track z0.
Definition: TTTrack.h:330
unsigned int tidRing(const DetId &id) const
Log< level::Warning, false > LogWarning
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511