CMS 3D CMS Logo

L1TrackJetProducer.cc
Go to the documentation of this file.
1 // Original simulation Author: Rishi Patel
2 //
3 // Rewritting/improvements: George Karathanasis,
4 // georgios.karathanasis@cern.ch, CU Boulder
5 //
6 // Created: Wed, 01 Aug 2018 14:01:41 GMT
7 // Latest update: Nov 2022 (by GK)
8 //
9 // Track jets are clustered in a two-layer process, first by clustering in phi,
10 // then by clustering in eta. The code proceeds as following: putting all tracks
11 // 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 // system include files
26 
27 // user include files
36 
37 //own headers
38 #include "L1TrackJetClustering.h"
39 
40 using namespace std;
41 using namespace edm;
42 using namespace l1t;
43 using namespace l1ttrackjet;
44 
46 public:
47  explicit L1TrackJetProducer(const ParameterSet &);
48  ~L1TrackJetProducer() override = default;
50  typedef vector<L1TTTrackType> L1TTTrackCollectionType;
52  static void fillDescriptions(ConfigurationDescriptions &descriptions);
53 
54 private:
55  void produce(Event &, const EventSetup &) override;
56 
57  // ----------member data ---------------------------
58 
59  vector<Ptr<L1TTTrackType>> L1TrkPtrs_;
60  vector<int> tdtrk_;
61  const float trkZMax_;
62  const float trkPtMax_;
63  const float trkPtMin_;
64  const float trkEtaMax_;
65  const float nStubs4PromptChi2_;
66  const float nStubs5PromptChi2_;
67  const float nStubs4PromptBend_;
68  const float nStubs5PromptBend_;
69  const int trkNPSStubMin_;
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  float zStep_;
79  float etaStep_;
80  float phiStep_;
81  const bool displaced_;
82  const float d0CutNStubs4_;
83  const float d0CutNStubs5_;
84  const float nStubs4DisplacedChi2_;
85  const float nStubs5DisplacedChi2_;
86  const float nStubs4DisplacedBend_;
87  const float nStubs5DisplacedBend_;
88  const int nDisplacedTracks_;
89  const float dzPVTrk_;
90 
94 };
95 
97  : trkZMax_(iConfig.getParameter<double>("trk_zMax")),
98  trkPtMax_(iConfig.getParameter<double>("trk_ptMax")),
99  trkPtMin_(iConfig.getParameter<double>("trk_ptMin")),
100  trkEtaMax_(iConfig.getParameter<double>("trk_etaMax")),
101  nStubs4PromptChi2_(iConfig.getParameter<double>("nStubs4PromptChi2")),
102  nStubs5PromptChi2_(iConfig.getParameter<double>("nStubs5PromptChi2")),
103  nStubs4PromptBend_(iConfig.getParameter<double>("nStubs4PromptBend")),
104  nStubs5PromptBend_(iConfig.getParameter<double>("nStubs5PromptBend")),
105  trkNPSStubMin_(iConfig.getParameter<int>("trk_nPSStubMin")),
106  lowpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
107  lowpTJetThreshold_(iConfig.getParameter<double>("lowpTJetThreshold")),
108  highpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
109  highpTJetThreshold_(iConfig.getParameter<double>("highpTJetThreshold")),
110  zBins_(iConfig.getParameter<int>("zBins")),
111  etaBins_(iConfig.getParameter<int>("etaBins")),
112  phiBins_(iConfig.getParameter<int>("phiBins")),
113  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
114  displaced_(iConfig.getParameter<bool>("displaced")),
115  d0CutNStubs4_(iConfig.getParameter<double>("d0_cutNStubs4")),
116  d0CutNStubs5_(iConfig.getParameter<double>("d0_cutNStubs5")),
117  nStubs4DisplacedChi2_(iConfig.getParameter<double>("nStubs4DisplacedChi2")),
118  nStubs5DisplacedChi2_(iConfig.getParameter<double>("nStubs5DisplacedChi2")),
119  nStubs4DisplacedBend_(iConfig.getParameter<double>("nStubs4DisplacedBend")),
120  nStubs5DisplacedBend_(iConfig.getParameter<double>("nStubs5DisplacedBend")),
121  nDisplacedTracks_(iConfig.getParameter<int>("nDisplacedTracks")),
122  dzPVTrk_(iConfig.getParameter<double>("MaxDzTrackPV")),
124  trackToken_(consumes<L1TTTrackRefCollectionType>(iConfig.getParameter<InputTag>("L1TrackInputTag"))),
125  PVtxToken_(consumes<l1t::VertexWordCollection>(iConfig.getParameter<InputTag>("L1PVertexInputTag"))) {
126  zStep_ = 2.0 * trkZMax_ / (zBins_ + 1); // added +1 in denom
127  etaStep_ = 2.0 * trkEtaMax_ / etaBins_; //etaStep is the width of an etabin
128  phiStep_ = 2 * M_PI / phiBins_;
129 
130  if (displaced_)
131  produces<TkJetCollection>("L1TrackJetsExtended");
132  else
133  produces<TkJetCollection>("L1TrackJets");
134 }
135 
137  unique_ptr<TkJetCollection> L1L1TrackJetProducer(new TkJetCollection);
138 
139  // Read inputs
140  const TrackerTopology &tTopo = iSetup.getData(tTopoToken_);
141 
143  iEvent.getByToken(trackToken_, TTTrackHandle);
144 
146  iEvent.getByToken(PVtxToken_, PVtx);
147  float PVz = (PVtx->at(0)).z0();
148 
149  L1TrkPtrs_.clear();
150  tdtrk_.clear();
151 
152  // track selection
153  for (unsigned int this_l1track = 0; this_l1track < TTTrackHandle->size(); this_l1track++) {
154  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
155  float trk_pt = trkPtr->momentum().perp();
156  int trk_nstubs = (int)trkPtr->getStubRefs().size();
157  float trk_chi2dof = trkPtr->chi2Red();
158  float trk_d0 = trkPtr->d0();
159  float trk_bendchi2 = trkPtr->stubPtConsistency();
160 
161  int trk_nPS = 0;
162  for (int istub = 0; istub < trk_nstubs; istub++) { // loop over the stubs
163  DetId detId(trkPtr->getStubRefs().at(istub)->getDetId());
164  if (detId.det() == DetId::Detector::Tracker) {
165  if ((detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3) ||
166  (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9))
167  trk_nPS++;
168  }
169  }
170  // select tracks
171  if (trk_nPS < trkNPSStubMin_)
172  continue;
173  if (!TrackQualitySelection(trk_nstubs,
174  trk_chi2dof,
175  trk_bendchi2,
184  displaced_))
185  continue;
186  if (std::abs(PVz - trkPtr->z0()) > dzPVTrk_ && dzPVTrk_ > 0)
187  continue;
188  if (std::abs(trkPtr->z0()) > trkZMax_)
189  continue;
190  if (std::abs(trkPtr->momentum().eta()) > trkEtaMax_)
191  continue;
192  if (trk_pt < trkPtMin_)
193  continue;
194  L1TrkPtrs_.push_back(trkPtr);
195 
196  if ((std::abs(trk_d0) > d0CutNStubs5_ && trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
197  (trk_nstubs == 4 && std::abs(trk_d0) > d0CutNStubs4_ && d0CutNStubs4_ >= 0))
198  tdtrk_.push_back(1); //displaced track
199  else
200  tdtrk_.push_back(0); // not displaced track
201  }
202 
203  // if no tracks pass selection return empty containers
204  if (L1TrkPtrs_.empty()) {
205  if (displaced_)
206  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
207  else
208  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
209  return;
210  }
211 
212  MaxZBin mzb;
213  mzb.znum = -1;
214  mzb.nclust = -1;
215  mzb.ht = -1;
216 
217  // create 2D grid of eta phi bins
218  EtaPhiBin epbins_default[phiBins_][etaBins_];
219  float phi = -1.0 * M_PI;
220  for (int i = 0; i < phiBins_; ++i) {
221  float eta = -1.0 * trkEtaMax_;
222  for (int j = 0; j < etaBins_; ++j) {
223  epbins_default[i][j].phi = (phi + (phi + phiStep_)) / 2.0;
224  epbins_default[i][j].eta = (eta + (eta + etaStep_)) / 2.0;
225  eta += etaStep_;
226  } // for each etabin
227  phi += phiStep_;
228  } // for each phibin (finished creating bins)
229 
230  // create z-bins (might be useful for displaced if we run w/o dz cut)
231  std::vector<float> zmins, zmaxs;
232  for (int zbin = 0; zbin < zBins_; zbin++) {
233  zmins.push_back(-1.0 * trkZMax_ + zStep_ * zbin);
234  zmaxs.push_back(-1.0 * trkZMax_ + zStep_ * zbin + zStep_ * 2);
235  }
236 
237  // create vectors of clusters
238  std::vector<std::vector<EtaPhiBin>> L1clusters;
239  L1clusters.reserve(phiBins_);
240  std::vector<EtaPhiBin> L2clusters;
241 
242  // default (empty) grid
243  for (int i = 0; i < phiBins_; ++i) {
244  for (int j = 0; j < etaBins_; ++j) {
245  epbins_default[i][j].pTtot = 0;
246  epbins_default[i][j].used = false;
247  epbins_default[i][j].ntracks = 0;
248  epbins_default[i][j].nxtracks = 0;
249  epbins_default[i][j].trackidx.clear();
250  }
251  }
252 
253  for (unsigned int zbin = 0; zbin < zmins.size(); ++zbin) {
254  // initialize grid
255  float zmin = zmins[zbin];
256  float zmax = zmaxs[zbin];
257  EtaPhiBin epbins[phiBins_][etaBins_];
258  std::copy(&epbins_default[0][0], &epbins_default[0][0] + phiBins_ * etaBins_, &epbins[0][0]);
259 
260  //clear cluster containers
261  L1clusters.clear();
262  L2clusters.clear();
263 
264  // fill grid with tracks
265  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
266  float trkZ = L1TrkPtrs_[k]->z0();
267  if (zmax < trkZ)
268  continue;
269  if (zmin > trkZ)
270  continue;
271  if (zbin == 0 && zmin == trkZ)
272  continue;
273  float trkpt = L1TrkPtrs_[k]->momentum().perp();
274  float trketa = L1TrkPtrs_[k]->momentum().eta();
275  float trkphi = L1TrkPtrs_[k]->momentum().phi();
276  for (int i = 0; i < phiBins_; ++i) {
277  for (int j = 0; j < etaBins_; ++j) {
278  float eta_min = epbins[i][j].eta - etaStep_ / 2.0; //eta min
279  float eta_max = epbins[i][j].eta + etaStep_ / 2.0; //eta max
280  float phi_min = epbins[i][j].phi - phiStep_ / 2.0; //phi min
281  float phi_max = epbins[i][j].phi + phiStep_ / 2.0; //phi max
282 
283  if ((trketa < eta_min) || (trketa > eta_max) || (trkphi < phi_min) || (trkphi > phi_max))
284  continue;
285 
286  if (trkpt < trkPtMax_)
287  epbins[i][j].pTtot += trkpt;
288  else
289  epbins[i][j].pTtot += trkPtMax_;
290  epbins[i][j].nxtracks += tdtrk_[k];
291  epbins[i][j].trackidx.push_back(k);
292  ++epbins[i][j].ntracks;
293  } // for each etabin
294  } // for each phibin
295  } //end loop over tracks
296 
297  // cluster tracks in eta (first layer) using grid
298  for (int phibin = 0; phibin < phiBins_; ++phibin) {
299  L1clusters.push_back(L1_clustering<EtaPhiBin, float, float, float>(epbins[phibin], etaBins_, etaStep_));
300  }
301 
302  // second layer clustering in phi for using eta clusters
303  L2clusters = L2_clustering<EtaPhiBin, float, float, float>(L1clusters, phiBins_, phiStep_, etaStep_);
304 
305  // sum all cluster pTs in this zbin to find max
306  float sum_pt = 0;
307  for (unsigned int k = 0; k < L2clusters.size(); ++k) {
308  if (L2clusters[k].pTtot > lowpTJetThreshold_ && L2clusters[k].ntracks < lowpTJetMinTrackMultiplicity_)
309  continue;
310  if (L2clusters[k].pTtot > highpTJetThreshold_ && L2clusters[k].ntracks < highpTJetMinTrackMultiplicity_)
311  continue;
312  if (L2clusters[k].pTtot > minTrkJetpT_)
313  sum_pt += L2clusters[k].pTtot;
314  }
315 
316  if (sum_pt < mzb.ht)
317  continue;
318 
319  // if ht is larger than previous max, this is the new vertex zbin
320  mzb.ht = sum_pt;
321  mzb.znum = zbin;
322  mzb.clusters = L2clusters;
323  mzb.nclust = L2clusters.size();
324  mzb.zbincenter = (zmin + zmax) / 2.0;
325  } //zbin loop
326 
327  // output
328  vector<Ptr<L1TTTrackType>> L1TrackAssocJet;
329  for (unsigned int j = 0; j < mzb.clusters.size(); ++j) {
330  float jetEta = mzb.clusters[j].eta;
331  float jetPhi = mzb.clusters[j].phi;
332  float jetPt = mzb.clusters[j].pTtot;
333  float jetPx = jetPt * cos(jetPhi);
334  float jetPy = jetPt * sin(jetPhi);
335  float jetPz = jetPt * sinh(jetEta);
336  float jetP = jetPt * cosh(jetEta);
337  int totalDisptrk = mzb.clusters[j].nxtracks;
338  bool isDispJet = (totalDisptrk > nDisplacedTracks_ || totalDisptrk == nDisplacedTracks_);
339 
340  math::XYZTLorentzVector jetP4(jetPx, jetPy, jetPz, jetP);
341  L1TrackAssocJet.clear();
342  for (unsigned int itrk = 0; itrk < mzb.clusters[j].trackidx.size(); itrk++)
343  L1TrackAssocJet.push_back(L1TrkPtrs_[mzb.clusters[j].trackidx[itrk]]);
344 
345  TkJet trkJet(jetP4, L1TrackAssocJet, mzb.zbincenter, mzb.clusters[j].ntracks, 0, totalDisptrk, 0, isDispJet);
346 
347  L1L1TrackJetProducer->push_back(trkJet);
348  }
349 
350  std::sort(
351  L1L1TrackJetProducer->begin(), L1L1TrackJetProducer->end(), [](auto &a, auto &b) { return a.pt() > b.pt(); });
352  if (displaced_)
353  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
354  else
355  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
356 }
357 
360  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("l1tTTTracksFromTrackletEmulation", "Level1TTTracks"));
361  desc.add<edm::InputTag>("L1PVertexInputTag", edm::InputTag("l1tVertexFinderEmulator", "L1VerticesEmulation"));
362  desc.add<double>("MaxDzTrackPV", 1.0);
363  desc.add<double>("trk_zMax", 15.0);
364  desc.add<double>("trk_ptMax", 200.0);
365  desc.add<double>("trk_ptMin", 3.0);
366  desc.add<double>("trk_etaMax", 2.4);
367  desc.add<double>("nStubs4PromptChi2", 5.0);
368  desc.add<double>("nStubs4PromptBend", 1.7);
369  desc.add<double>("nStubs5PromptChi2", 2.75);
370  desc.add<double>("nStubs5PromptBend", 3.5);
371  desc.add<int>("trk_nPSStubMin", -1);
372  desc.add<double>("minTrkJetpT", -1.0);
373  desc.add<int>("etaBins", 24);
374  desc.add<int>("phiBins", 27);
375  desc.add<int>("zBins", 1);
376  desc.add<double>("d0_cutNStubs4", -1);
377  desc.add<double>("d0_cutNStubs5", -1);
378  desc.add<int>("lowpTJetMinTrackMultiplicity", 2);
379  desc.add<double>("lowpTJetThreshold", 50.0);
380  desc.add<int>("highpTJetMinTrackMultiplicity", 3);
381  desc.add<double>("highpTJetThreshold", 100.0);
382  desc.add<bool>("displaced", false);
383  desc.add<double>("nStubs4DisplacedChi2", 5.0);
384  desc.add<double>("nStubs4DisplacedBend", 1.7);
385  desc.add<double>("nStubs5DisplacedChi2", 2.75);
386  desc.add<double>("nStubs5DisplacedBend", 3.5);
387  desc.add<int>("nDisplacedTracks", 2);
388  descriptions.add("l1tTrackJets", desc);
389 }
390 
391 //define this as a plug-in
const int highpTJetMinTrackMultiplicity_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int tobLayer(const DetId &id) const
const float lowpTJetThreshold_
T perp() const
Definition: PV3DBase.h:69
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
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_)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< unsigned int > trackidx
T eta() const
Definition: PV3DBase.h:73
const float nStubs4DisplacedBend_
delete x;
Definition: CaloConfig.h:22
L1TrackJetProducer(const ParameterSet &)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const float nStubs4DisplacedChi2_
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
const float nStubs5PromptChi2_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const float nStubs4PromptChi2_
int iEvent
Definition: GenABIO.cc:224
vector< L1TTTrackType > L1TTTrackCollectionType
static void fillDescriptions(ConfigurationDescriptions &descriptions)
const float highpTJetThreshold_
const float nStubs5DisplacedBend_
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
std::vector< EtaPhiBin > clusters
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const float nStubs4PromptBend_
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:197
static constexpr auto TOB
const float nStubs5PromptBend_
const float nStubs5DisplacedChi2_
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
const EDGetTokenT< l1t::VertexWordCollection > PVtxToken_
Definition: DetId.h:17
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
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 b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double d0() const
Track d0.
Definition: TTTrack.h:325
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
const int lowpTJetMinTrackMultiplicity_
HLT enums.
double chi2Red() const
Chi2 reduced.
Definition: TTTrack.h:359
double a
Definition: hdecay.h:121
double z0() const
Track z0.
Definition: TTTrack.h:330
unsigned int tidRing(const DetId &id) const
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
edm::RefVector< L1TTTrackCollectionType > L1TTTrackRefCollectionType
def move(src, dest)
Definition: eostools.py:511