CMS 3D CMS Logo

Phase2L1CaloJetEmulator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/L1CaloTrigger
4 // Class: Phase2L1CaloJetEmulator
5 //
13 //
14 // Original Author: Pallabi Das
15 // Created: Tue, 11 Apr 2023 11:27:33 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
31 
43 
45 #include <ap_int.h>
46 #include <cstdio>
47 #include <fstream>
48 #include <iomanip>
49 #include <iostream>
50 #include "TF1.h"
51 
52 //
53 // class declaration
54 //
55 
57 public:
59  ~Phase2L1CaloJetEmulator() override;
60 
61  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
62 
63 private:
64  void produce(edm::Event&, const edm::EventSetup&) override;
65  float get_jet_pt_calibration(const float& jet_pt, const float& jet_eta) const;
66  float get_tau_pt_calibration(const float& tau_pt, const float& tau_eta) const;
67 
68  // ----------member data ---------------------------
73  std::vector<edm::ParameterSet> nHits_to_nvtx_params;
74  std::vector<edm::ParameterSet> nvtx_to_PU_sub_params;
75  std::map<std::string, TF1> nHits_to_nvtx_funcs;
76  std::map<std::string, TF1> hgcalEM_nvtx_to_PU_sub_funcs;
77  std::map<std::string, TF1> hgcalHad_nvtx_to_PU_sub_funcs;
78  std::map<std::string, TF1> hf_nvtx_to_PU_sub_funcs;
79  std::map<std::string, std::map<std::string, TF1>> all_nvtx_to_PU_sub_funcs;
80 
81  // For fetching jet pt calibrations
82  std::vector<double> jetPtBins;
83  std::vector<double> absEtaBinsBarrel;
84  std::vector<double> jetCalibrationsBarrel;
85  std::vector<double> absEtaBinsHGCal;
86  std::vector<double> jetCalibrationsHGCal;
87  std::vector<double> absEtaBinsHF;
88  std::vector<double> jetCalibrationsHF;
89 
90  // For fetching tau pt calibrations
91  std::vector<double> tauPtBins;
92  std::vector<double> tauAbsEtaBinsBarrel;
93  std::vector<double> tauCalibrationsBarrel;
94  std::vector<double> tauAbsEtaBinsHGCal;
95  std::vector<double> tauCalibrationsHGCal;
96 
97  // For storing jet calibrations
98  std::vector<std::vector<double>> calibrationsBarrel;
99  std::vector<std::vector<double>> calibrationsHGCal;
100  std::vector<std::vector<double>> calibrationsHF;
101 
102  // For storing tau calibrations
103  std::vector<std::vector<double>> tauPtCalibrationsBarrel;
104  std::vector<std::vector<double>> tauPtCalibrationsHGCal;
105 };
106 
107 //
108 // constructors and destructor
109 //
111  : caloTowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("gctFullTowers"))),
112  hgcalTowerToken_(consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("hgcalTowers"))),
113  hfToken_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("hcalDigis"))),
114  decoderTag_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>(edm::ESInputTag("", ""))),
115  nHits_to_nvtx_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nHits_to_nvtx_params")),
116  nvtx_to_PU_sub_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nvtx_to_PU_sub_params")),
117  jetPtBins(iConfig.getParameter<std::vector<double>>("jetPtBins")),
118  absEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("absEtaBinsBarrel")),
119  jetCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("jetCalibrationsBarrel")),
120  absEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("absEtaBinsHGCal")),
121  jetCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("jetCalibrationsHGCal")),
122  absEtaBinsHF(iConfig.getParameter<std::vector<double>>("absEtaBinsHF")),
123  jetCalibrationsHF(iConfig.getParameter<std::vector<double>>("jetCalibrationsHF")),
124  tauPtBins(iConfig.getParameter<std::vector<double>>("tauPtBins")),
125  tauAbsEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsBarrel")),
126  tauCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("tauCalibrationsBarrel")),
127  tauAbsEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsHGCal")),
128  tauCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("tauCalibrationsHGCal")) {
129  for (uint i = 0; i < nHits_to_nvtx_params.size(); i++) {
131  std::string calo = pset->getParameter<std::string>("fit");
132  nHits_to_nvtx_funcs[calo.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
133  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(0, pset->getParameter<std::vector<double>>("nHits_params").at(0));
134  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(1, pset->getParameter<std::vector<double>>("nHits_params").at(1));
135  }
139  for (uint i = 0; i < nvtx_to_PU_sub_params.size(); i++) {
141  std::string calo = pset->getParameter<std::string>("calo");
142  std::string iEta = pset->getParameter<std::string>("iEta");
143  double p1 = pset->getParameter<std::vector<double>>("nvtx_params").at(0);
144  double p2 = pset->getParameter<std::vector<double>>("nvtx_params").at(1);
145 
146  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
147  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(0, p1);
148  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(1, p2);
149  }
150 
151  // Fill the jet pt calibration 2D vector
152  // Dimension 1 is AbsEta bin
153  // Dimension 2 is jet pT bin which is filled with the actual callibration value
154  // size()-1 b/c the inputs have lower and upper bounds
155  // Do Barrel, then HGCal, then HF
156  int index = 0;
157  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsBarrel.size() - 1; abs_eta++) {
158  std::vector<double> pt_bin_calibs;
159  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
160  pt_bin_calibs.push_back(jetCalibrationsBarrel.at(index));
161  index++;
162  }
163  calibrationsBarrel.push_back(pt_bin_calibs);
164  }
165 
166  index = 0;
167  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHGCal.size() - 1; abs_eta++) {
168  std::vector<double> pt_bin_calibs;
169  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
170  pt_bin_calibs.push_back(jetCalibrationsHGCal.at(index));
171  index++;
172  }
173  calibrationsHGCal.push_back(pt_bin_calibs);
174  }
175 
176  index = 0;
177  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHF.size() - 1; abs_eta++) {
178  std::vector<double> pt_bin_calibs;
179  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
180  pt_bin_calibs.push_back(jetCalibrationsHF.at(index));
181  index++;
182  }
183  calibrationsHF.push_back(pt_bin_calibs);
184  }
185 
186  // Fill the tau pt calibration 2D vector
187  // Dimension 1 is AbsEta bin
188  // Dimension 2 is tau pT bin which is filled with the actual calibration value
189  // Do Barrel, then HGCal
190  //
191  // Note to future developers: be very concious of the order in which the calibrations are printed
192  // out in tool which makse the cfg files. You need to match that exactly when loading them and
193  // using the calibrations below.
194  index = 0;
195  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsBarrel.size() - 1; abs_eta++) {
196  std::vector<double> pt_bin_calibs;
197  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
198  pt_bin_calibs.push_back(tauCalibrationsBarrel.at(index));
199  index++;
200  }
201  tauPtCalibrationsBarrel.push_back(pt_bin_calibs);
202  }
203 
204  index = 0;
205  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsHGCal.size() - 1; abs_eta++) {
206  std::vector<double> pt_bin_calibs;
207  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
208  pt_bin_calibs.push_back(tauCalibrationsHGCal.at(index));
209  index++;
210  }
211  tauPtCalibrationsHGCal.push_back(pt_bin_calibs);
212  }
213 
214  produces<l1tp2::Phase2L1CaloJetCollection>("GCTJet");
215 }
216 
218 
219 //
220 // member functions
221 //
222 
223 // ------------ method called to produce the data ------------
225  using namespace edm;
226  std::unique_ptr<l1tp2::Phase2L1CaloJetCollection> jetCands(make_unique<l1tp2::Phase2L1CaloJetCollection>());
227 
228  // Assign ETs to each eta-half of the barrel region (17x72 --> 18x72 to be able to make 3x3 super towers)
229  edm::Handle<std::vector<l1tp2::CaloTower>> caloTowerCollection;
230  if (!iEvent.getByToken(caloTowerToken_, caloTowerCollection))
231  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from caloTowerCollection!";
232 
233  iEvent.getByToken(caloTowerToken_, caloTowerCollection);
234  float GCTintTowers[nBarrelEta][nBarrelPhi];
235  float realEta[nBarrelEta][nBarrelPhi];
236  float realPhi[nBarrelEta][nBarrelPhi];
237  for (const l1tp2::CaloTower& i : *caloTowerCollection) {
238  int ieta = i.towerIEta();
239  int iphi = i.towerIPhi();
240  if (i.ecalTowerEt() > 1.)
241  GCTintTowers[ieta][iphi] = i.ecalTowerEt(); // suppress <= 1 GeV towers
242  else
243  GCTintTowers[ieta][iphi] = 0;
244  realEta[ieta][iphi] = i.towerEta();
245  realPhi[ieta][iphi] = i.towerPhi();
246  }
247 
248  float temporary[nBarrelEta / 2][nBarrelPhi];
249 
250  // HGCal and HF info used for nvtx estimation
251  edm::Handle<l1t::HGCalTowerBxCollection> hgcalTowerCollection;
252  if (!iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection))
253  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from hgcalTowerCollection!";
254  l1t::HGCalTowerBxCollection hgcalTowerColl;
255  iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection);
256  hgcalTowerColl = (*hgcalTowerCollection.product());
257 
259  if (!iEvent.getByToken(hfToken_, hfHandle))
260  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get HcalTrigPrimDigi for HF!";
261  iEvent.getByToken(hfToken_, hfHandle);
262 
263  int i_hgcalEM_hits_leq_threshold = 0;
264  int i_hgcalHad_hits_leq_threshold = 0;
265  int i_hf_hits_leq_threshold = 0;
266  for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
267  if (it->etEm() <= 1.75 && it->etEm() >= 1.25) {
268  i_hgcalEM_hits_leq_threshold++;
269  }
270  if (it->etHad() <= 1.25 && it->etHad() >= 0.75) {
271  i_hgcalHad_hits_leq_threshold++;
272  }
273  }
274  const auto& decoder = iSetup.getData(decoderTag_);
275  for (const auto& hit : *hfHandle.product()) {
276  double et = decoder.hcaletValue(hit.id(), hit.t0());
277  if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
278  continue;
279  if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
280  continue;
281  if (et <= 15.0 && et >= 10.0)
282  i_hf_hits_leq_threshold++;
283  }
284 
285  double hgcalEM_nvtx = nHits_to_nvtx_funcs["hgcalEM"].Eval(i_hgcalEM_hits_leq_threshold);
286  if (hgcalEM_nvtx < 0)
287  hgcalEM_nvtx = 0;
288  double hgcalHad_nvtx = nHits_to_nvtx_funcs["hgcalHad"].Eval(i_hgcalHad_hits_leq_threshold);
289  if (hgcalHad_nvtx < 0)
290  hgcalHad_nvtx = 0;
291  double hf_nvtx = nHits_to_nvtx_funcs["hf"].Eval(i_hf_hits_leq_threshold);
292  if (hf_nvtx < 0)
293  hf_nvtx = 0;
294  double EstimatedNvtx = (hgcalEM_nvtx + hgcalHad_nvtx + hf_nvtx) / 3.;
295 
296  // Assign ETs to each eta-half of the endcap region (18x72)
298  float hgcalEta[nHgcalEta][nHgcalPhi];
299  float hgcalPhi[nHgcalEta][nHgcalPhi];
300 
301  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
302  for (int ieta = 0; ieta < nHgcalEta; ieta++) {
303  hgcalTowers[ieta][iphi] = 0;
304  if (ieta < nHgcalEta / 2)
305  hgcalEta[ieta][iphi] = -3.045 + ieta * 0.087 + 0.0435;
306  else
307  hgcalEta[ieta][iphi] = 1.479 + (ieta - nHgcalEta / 2) * 0.087 + 0.0435;
308  hgcalPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHgcalPhi) + (M_PI / nHgcalPhi);
309  }
310  }
311 
312  for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
313  float eta = it->eta();
314  int ieta;
315  if (eta < 0)
316  ieta = 19 - it->id().iEta();
317  else
318  ieta = 20 + it->id().iEta();
319  if (eta > 1.479)
320  ieta = ieta - 4;
321  int iphi = it->id().iPhi();
322 
323  float hgcal_etEm = it->etEm();
324  float hgcal_etHad = it->etHad();
325  std::string etaKey = "";
326  if (abs(eta) <= 1.8)
327  etaKey = "er1p4to1p8";
328  else if (abs(eta) <= 2.1 && abs(eta) > 1.8)
329  etaKey = "er1p8to2p1";
330  else if (abs(eta) <= 2.4 && abs(eta) > 2.1)
331  etaKey = "er2p1to2p4";
332  else if (abs(eta) <= 2.7 && abs(eta) > 2.4)
333  etaKey = "er2p4to2p7";
334  else if (abs(eta) <= 3.1 && abs(eta) > 2.7)
335  etaKey = "er2p7to3p1";
336  if (!etaKey.empty()) {
337  hgcal_etEm = it->etEm() - all_nvtx_to_PU_sub_funcs["hgcalEM"][etaKey].Eval(EstimatedNvtx);
338  hgcal_etHad = it->etHad() - all_nvtx_to_PU_sub_funcs["hgcalHad"][etaKey].Eval(EstimatedNvtx);
339  }
340 
341  if (hgcal_etEm < 0)
342  hgcal_etEm = 0;
343  if (hgcal_etHad < 0)
344  hgcal_etHad = 0;
345  if ((it->etEm() + it->etHad() > 1.) && abs(eta) > 1.479)
346  hgcalTowers[ieta][iphi] = hgcal_etEm + hgcal_etHad; // suppress <= 1 GeV towers
347  }
348 
349  float temporary_hgcal[nHgcalEta / 2][nHgcalPhi];
350 
351  // Assign ETs to each eta-half of the forward region (12x72)
352  float hfTowers[nHfEta][nHfPhi];
353  float hfEta[nHfEta][nHfPhi];
354  float hfPhi[nHfEta][nHfPhi];
355  for (int iphi = 0; iphi < nHfPhi; iphi++) {
356  for (int ieta = 0; ieta < nHfEta; ieta++) {
357  hfTowers[ieta][iphi] = 0;
358  int temp;
359  if (ieta < nHfEta / 2)
361  else
364  hfPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHfPhi) + (M_PI / nHfPhi);
365  }
366  }
367 
368  for (const auto& hit : *hfHandle.product()) {
369  double et = decoder.hcaletValue(hit.id(), hit.t0());
370  int ieta = 0;
371  if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
372  continue;
373  if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
374  continue;
375  if (hit.id().ieta() <= -(l1t::CaloTools::kHFBegin + 1)) {
376  ieta = hit.id().ieta() + l1t::CaloTools::kHFEnd;
377  } else if (hit.id().ieta() >= (l1t::CaloTools::kHFBegin + 1)) {
378  ieta = nHfEta / 2 + (hit.id().ieta() - (l1t::CaloTools::kHFBegin + 1));
379  }
380  int iphi = 0;
381  if (hit.id().iphi() <= nHfPhi / 2)
382  iphi = hit.id().iphi() + (nHfPhi / 2 - 1);
383  else if (hit.id().iphi() > nHfPhi / 2)
384  iphi = hit.id().iphi() - (nHfPhi / 2 + 1);
385  if (abs(hit.id().ieta()) <= 33 && abs(hit.id().ieta()) >= 29)
386  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er29to33"].Eval(EstimatedNvtx);
387  if (abs(hit.id().ieta()) <= 37 && abs(hit.id().ieta()) >= 34)
388  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er34to37"].Eval(EstimatedNvtx);
389  if (abs(hit.id().ieta()) <= 41 && abs(hit.id().ieta()) >= 38)
390  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er38to41"].Eval(EstimatedNvtx);
391  if (et < 0)
392  et = 0;
393  if (et > 1.)
394  hfTowers[ieta][iphi] = et; // suppress <= 1 GeV towers
395  }
396 
397  float temporary_hf[nHfEta / 2][nHfPhi];
398 
399  // Begin creating jets
400  // First create 3x3 super towers: 6x24 in barrel, endcap; 4x24 in forward
401  // Then create up to 10 jets in each eta half of barrel, endcap, forward regions
402 
403  vector<l1tp2::Phase2L1CaloJet> halfBarrelJets, halfHgcalJets, halfHfJets;
404  halfBarrelJets.clear();
405  halfHgcalJets.clear();
406  halfHfJets.clear();
407  vector<l1tp2::Phase2L1CaloJet> allJets;
408  allJets.clear();
409 
410  for (int k = 0; k < 2; k++) {
411  halfBarrelJets.clear();
412  halfHgcalJets.clear();
413  halfHfJets.clear();
415 
416  // BARREL
417  for (int iphi = 0; iphi < nBarrelPhi; iphi++) {
418  for (int ieta = 0; ieta < nBarrelEta / 2; ieta++) {
419  if (k == 0)
420  temporary[ieta][iphi] = GCTintTowers[ieta][iphi];
421  else
422  temporary[ieta][iphi] = GCTintTowers[nBarrelEta / 2 + ieta][iphi];
423  }
424  }
425 
427  gctobj::makeST(temporary, tempST);
428  float TTseedThresholdBarrel = 5.;
429 
430  for (int i = 0; i < nJets; i++) {
431  jet[i] = gctobj::getRegion(tempST, TTseedThresholdBarrel);
432  l1tp2::Phase2L1CaloJet tempJet;
433  int gctjeteta = jet[i].etaCenter;
434  int gctjetphi = jet[i].phiCenter;
435  tempJet.setJetIEta(gctjeteta + k * nBarrelEta / 2);
436  tempJet.setJetIPhi(gctjetphi);
437  float jeteta = realEta[gctjeteta + k * nBarrelEta / 2][gctjetphi];
438  float jetphi = realPhi[gctjeteta + k * nBarrelEta / 2][gctjetphi];
439  tempJet.setJetEta(jeteta);
440  tempJet.setJetPhi(jetphi);
441  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
442  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
443  tempJet.setTowerEt(jet[i].energyMax);
444  int gcttowereta = jet[i].etaMax;
445  int gcttowerphi = jet[i].phiMax;
446  tempJet.setTowerIEta(gcttowereta + k * nBarrelEta / 2);
447  tempJet.setTowerIPhi(gcttowerphi);
448  float towereta = realEta[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
449  float towerphi = realPhi[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
450  tempJet.setTowerEta(towereta);
451  tempJet.setTowerPhi(towerphi);
453  tempJetp4.SetPt(tempJet.jetEt());
454  tempJetp4.SetEta(tempJet.jetEta());
455  tempJetp4.SetPhi(tempJet.jetPhi());
456  tempJetp4.SetM(0.);
457  tempJet.setP4(tempJetp4);
458 
459  if (jet[i].energy > 0.)
460  halfBarrelJets.push_back(tempJet);
461  }
462 
463  // ENDCAP
464  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
465  for (int ieta = 0; ieta < nHgcalEta / 2; ieta++) {
466  if (k == 0)
467  temporary_hgcal[ieta][iphi] = hgcalTowers[ieta][iphi];
468  else
469  temporary_hgcal[ieta][iphi] = hgcalTowers[nHgcalEta / 2 + ieta][iphi];
470  }
471  }
472 
473  gctobj::GCTsupertower_t tempST_hgcal[nSTEta][nSTPhi];
474  gctobj::makeST_hgcal(temporary_hgcal, tempST_hgcal);
475  float TTseedThresholdEndcap = 3.;
476  for (int i = nJets; i < 2 * nJets; i++) {
477  jet[i] = gctobj::getRegion(tempST_hgcal, TTseedThresholdEndcap);
478  l1tp2::Phase2L1CaloJet tempJet;
479  int hgcaljeteta = jet[i].etaCenter;
480  int hgcaljetphi = jet[i].phiCenter;
481  tempJet.setJetIEta(hgcaljeteta + k * nHgcalEta / 2);
482  tempJet.setJetIPhi(hgcaljetphi);
483  float jeteta = hgcalEta[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
484  float jetphi = hgcalPhi[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
485  tempJet.setJetEta(jeteta);
486  tempJet.setJetPhi(jetphi);
487  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
488  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
489  tempJet.setTowerEt(jet[i].energyMax);
490  int hgcaltowereta = jet[i].etaMax;
491  int hgcaltowerphi = jet[i].phiMax;
492  tempJet.setTowerIEta(hgcaltowereta + k * nHgcalEta / 2);
493  tempJet.setTowerIPhi(hgcaltowerphi);
494  float towereta = hgcalEta[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
495  float towerphi = hgcalPhi[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
496  tempJet.setTowerEta(towereta);
497  tempJet.setTowerPhi(towerphi);
499  tempJetp4.SetPt(tempJet.jetEt());
500  tempJetp4.SetEta(tempJet.jetEta());
501  tempJetp4.SetPhi(tempJet.jetPhi());
502  tempJetp4.SetM(0.);
503  tempJet.setP4(tempJetp4);
504 
505  if (jet[i].energy > 0.)
506  halfHgcalJets.push_back(tempJet);
507  }
508 
509  // HF
510  for (int iphi = 0; iphi < nHfPhi; iphi++) {
511  for (int ieta = 0; ieta < nHfEta / 2; ieta++) {
512  if (k == 0)
513  temporary_hf[ieta][iphi] = hfTowers[ieta][iphi];
514  else
515  temporary_hf[ieta][iphi] = hfTowers[nHfEta / 2 + ieta][iphi];
516  }
517  }
518 
520  gctobj::makeST_hf(temporary_hf, tempST_hf);
521  float TTseedThresholdHF = 5.;
522  for (int i = 2 * nJets; i < 3 * nJets; i++) {
523  jet[i] = gctobj::getRegion(tempST_hf, TTseedThresholdHF);
524  l1tp2::Phase2L1CaloJet tempJet;
525  int hfjeteta = jet[i].etaCenter;
526  int hfjetphi = jet[i].phiCenter;
527  tempJet.setJetIEta(hfjeteta + k * nHfEta / 2);
528  tempJet.setJetIPhi(hfjetphi);
529  float jeteta = hfEta[hfjeteta + k * nHfEta / 2][hfjetphi];
530  float jetphi = hfPhi[hfjeteta + k * nHfEta / 2][hfjetphi];
531  tempJet.setJetEta(jeteta);
532  tempJet.setJetPhi(jetphi);
533  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
534  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
535  tempJet.setTowerEt(jet[i].energyMax);
536  int hftowereta = jet[i].etaMax;
537  int hftowerphi = jet[i].phiMax;
538  tempJet.setTowerIEta(hftowereta + k * nHfEta / 2);
539  tempJet.setTowerIPhi(hftowerphi);
540  float towereta = hfEta[hftowereta + k * nHfEta / 2][hftowerphi];
541  float towerphi = hfPhi[hftowereta + k * nHfEta / 2][hftowerphi];
542  tempJet.setTowerEta(towereta);
543  tempJet.setTowerPhi(towerphi);
545  tempJetp4.SetPt(tempJet.jetEt());
546  tempJetp4.SetEta(tempJet.jetEta());
547  tempJetp4.SetPhi(tempJet.jetPhi());
548  tempJetp4.SetM(0.);
549  tempJet.setP4(tempJetp4);
550 
551  if (jet[i].energy > 0.)
552  halfHfJets.push_back(tempJet);
553  }
554 
555  // Stitching:
556  // if the jet eta is at the boundary: for HB it should be within 0-1 in -ve eta, 32-33 in +ve eta; for HE it should be within 0-1/16-17 in -ve eta, 34-35/18-19 in +ve eta; for HF it should be within 10-11 in -ve eta, 12-13 in +ve eta
557  // then get the phi of that jet and check if there is a neighbouring jet with the same phi, then merge to the jet that has higher ET
558  // in both eta/phi allow a maximum of one tower between jet centers for stitching
559 
560  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
561  if (halfHgcalJets.at(i).jetIEta() >= (nHgcalEta / 2 - 2) && halfHgcalJets.at(i).jetIEta() < (nHgcalEta / 2 + 2)) {
562  float hgcal_ieta = k * nBarrelEta + halfHgcalJets.at(i).jetIEta();
563  for (size_t j = 0; j < halfBarrelJets.size(); j++) {
564  float barrel_ieta = nHgcalEta / 2 + halfBarrelJets.at(j).jetIEta();
565  if (abs(barrel_ieta - hgcal_ieta) <= 2 &&
566  abs(halfBarrelJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) <= 2) {
567  float totalet = halfBarrelJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
568  float totalTauEt = halfBarrelJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
569  if (halfBarrelJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
570  halfHgcalJets.at(i).setJetEt(0.);
571  halfHgcalJets.at(i).setTauEt(0.);
572  halfBarrelJets.at(j).setJetEt(totalet);
573  halfBarrelJets.at(j).setTauEt(totalTauEt);
575  tempJetp4.SetPt(totalet);
576  tempJetp4.SetEta(halfBarrelJets.at(j).jetEta());
577  tempJetp4.SetPhi(halfBarrelJets.at(j).jetPhi());
578  tempJetp4.SetM(0.);
579  halfBarrelJets.at(j).setP4(tempJetp4);
580  } else {
581  halfHgcalJets.at(i).setJetEt(totalet);
582  halfHgcalJets.at(i).setTauEt(totalTauEt);
583  halfBarrelJets.at(j).setJetEt(0.);
584  halfBarrelJets.at(j).setTauEt(0.);
586  tempJetp4.SetPt(totalet);
587  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
588  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
589  tempJetp4.SetM(0.);
590  halfHgcalJets.at(i).setP4(tempJetp4);
591  }
592  }
593  }
594  } else if (halfHgcalJets.at(i).jetIEta() < 2 || halfHgcalJets.at(i).jetIEta() >= (nHgcalEta - 2)) {
595  float hgcal_ieta = k * nBarrelEta + nHfEta / 2 + halfHgcalJets.at(i).jetIEta();
596  for (size_t j = 0; j < halfHfJets.size(); j++) {
597  float hf_ieta = k * nBarrelEta + k * nHgcalEta + halfHfJets.at(j).jetIEta();
598  if (abs(hgcal_ieta - hf_ieta) < 3 && abs(halfHfJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) < 3) {
599  float totalet = halfHfJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
600  float totalTauEt = halfHfJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
601  if (halfHfJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
602  halfHgcalJets.at(i).setJetEt(0.);
603  halfHgcalJets.at(i).setTauEt(0.);
604  halfHfJets.at(j).setJetEt(totalet);
605  halfHfJets.at(j).setTauEt(totalTauEt);
607  tempJetp4.SetPt(totalet);
608  tempJetp4.SetEta(halfHfJets.at(j).jetEta());
609  tempJetp4.SetPhi(halfHfJets.at(j).jetPhi());
610  tempJetp4.SetM(0.);
611  halfHfJets.at(j).setP4(tempJetp4);
612  } else {
613  halfHgcalJets.at(i).setJetEt(totalet);
614  halfHgcalJets.at(i).setTauEt(totalTauEt);
615  halfHfJets.at(j).setJetEt(0.);
616  halfHfJets.at(j).setTauEt(0.);
618  tempJetp4.SetPt(totalet);
619  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
620  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
621  tempJetp4.SetM(0.);
622  halfHgcalJets.at(i).setP4(tempJetp4);
623  }
624  }
625  }
626  }
627  }
628 
629  // Write 6 leading jets from each eta half
630 
631  vector<l1tp2::Phase2L1CaloJet> halfAllJets;
632  halfAllJets.clear();
633 
634  std::sort(halfBarrelJets.begin(), halfBarrelJets.end(), gctobj::compareByEt);
635  for (size_t i = 0; i < halfBarrelJets.size(); i++) {
636  if (halfBarrelJets.at(i).jetEt() > 0. && i < 6)
637  halfAllJets.push_back(halfBarrelJets.at(i));
638  }
639 
640  std::sort(halfHgcalJets.begin(), halfHgcalJets.end(), gctobj::compareByEt);
641  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
642  if (halfHgcalJets.at(i).jetEt() > 0. && i < 6)
643  halfAllJets.push_back(halfHgcalJets.at(i));
644  }
645 
646  std::sort(halfHfJets.begin(), halfHfJets.end(), gctobj::compareByEt);
647  for (size_t i = 0; i < halfHfJets.size(); i++) {
648  if (halfHfJets.at(i).jetEt() > 0. && i < 6)
649  halfAllJets.push_back(halfHfJets.at(i));
650  }
651 
652  std::sort(halfAllJets.begin(), halfAllJets.end(), gctobj::compareByEt);
653  for (size_t i = 0; i < halfAllJets.size(); i++) {
654  if (halfAllJets.at(i).jetEt() > 0. && i < 6)
655  allJets.push_back(halfAllJets.at(i));
656  }
657  }
658 
659  std::sort(allJets.begin(), allJets.end(), gctobj::compareByEt);
660  for (size_t i = 0; i < allJets.size(); i++) {
661  jetCands->push_back(allJets.at(i));
662  }
663 
664  iEvent.put(std::move(jetCands), "GCTJet");
665 }
666 
667 // Apply calibrations to HCAL energy based on Jet Eta, Jet pT
668 float Phase2L1CaloJetEmulator::get_jet_pt_calibration(const float& jet_pt, const float& jet_eta) const {
669  float abs_eta = std::abs(jet_eta);
670  float tmp_jet_pt = jet_pt;
671  if (tmp_jet_pt > 499)
672  tmp_jet_pt = 499;
673 
674  // Different indices sizes in different calo regions.
675  // Barrel...
676  size_t eta_index = 0;
677  size_t pt_index = 0;
678  float calib = 1.0;
679  if (abs_eta <= 1.5) {
680  // Start loop checking 2nd value
681  for (unsigned int i = 1; i < absEtaBinsBarrel.size(); i++) {
682  if (abs_eta <= absEtaBinsBarrel.at(i))
683  break;
684  eta_index++;
685  }
686  // Start loop checking 2nd value
687  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
688  if (tmp_jet_pt <= jetPtBins.at(i))
689  break;
690  pt_index++;
691  }
692  calib = calibrationsBarrel[eta_index][pt_index];
693  } // end Barrel
694  else if (abs_eta <= 3.0) // HGCal
695  {
696  // Start loop checking 2nd value
697  for (unsigned int i = 1; i < absEtaBinsHGCal.size(); i++) {
698  if (abs_eta <= absEtaBinsHGCal.at(i))
699  break;
700  eta_index++;
701  }
702  // Start loop checking 2nd value
703  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
704  if (tmp_jet_pt <= jetPtBins.at(i))
705  break;
706  pt_index++;
707  }
708  calib = calibrationsHGCal[eta_index][pt_index];
709  } // end HGCal
710  else // HF
711  {
712  // Start loop checking 2nd value
713  for (unsigned int i = 1; i < absEtaBinsHF.size(); i++) {
714  if (abs_eta <= absEtaBinsHF.at(i))
715  break;
716  eta_index++;
717  }
718  // Start loop checking 2nd value
719  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
720  if (tmp_jet_pt <= jetPtBins.at(i))
721  break;
722  pt_index++;
723  }
724  calib = calibrationsHF[eta_index][pt_index];
725  } // end HF
726 
727  return jet_pt * calib;
728 }
729 
730 // Apply calibrations to tau pT based on L1EG info, EM Fraction, Tau Eta, Tau pT
731 float Phase2L1CaloJetEmulator::get_tau_pt_calibration(const float& tau_pt, const float& tau_eta) const {
732  float abs_eta = std::abs(tau_eta);
733  float tmp_tau_pt = tau_pt;
734  if (tmp_tau_pt > 199)
735  tmp_tau_pt = 199;
736 
737  // Different indices sizes in different calo regions.
738  // Barrel...
739  size_t eta_index = 0;
740  size_t pt_index = 0;
741  float calib = 1.0;
742  if (abs_eta <= 1.5) {
743  // Start loop checking 2nd value
744  for (unsigned int i = 1; i < tauAbsEtaBinsBarrel.size(); i++) {
745  if (abs_eta <= tauAbsEtaBinsBarrel.at(i))
746  break;
747  eta_index++;
748  }
749  // Start loop checking 2nd value
750  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
751  if (tmp_tau_pt <= tauPtBins.at(i))
752  break;
753  pt_index++;
754  }
755  calib = tauPtCalibrationsBarrel[eta_index][pt_index];
756  } // end Barrel
757  else if (abs_eta <= 3.0) // HGCal
758  {
759  // Start loop checking 2nd value
760  for (unsigned int i = 1; i < tauAbsEtaBinsHGCal.size(); i++) {
761  if (abs_eta <= tauAbsEtaBinsHGCal.at(i))
762  break;
763  eta_index++;
764  }
765  // Start loop checking 2nd value
766  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
767  if (tmp_tau_pt <= tauPtBins.at(i))
768  break;
769  pt_index++;
770  }
771  calib = tauPtCalibrationsHGCal[eta_index][pt_index];
772  } // end HGCal
773 
774  return tau_pt * calib;
775 }
776 
777 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
780  desc.add<edm::InputTag>("gctFullTowers", edm::InputTag("l1tPhase2L1CaloEGammaEmulator", "GCTFullTowers"));
781  desc.add<edm::InputTag>("hgcalTowers", edm::InputTag("l1tHGCalTowerProducer", "HGCalTowerProcessor"));
782  desc.add<edm::InputTag>("hcalDigis", edm::InputTag("simHcalTriggerPrimitiveDigis"));
783 
784  edm::ParameterSetDescription nHits_params_validator;
785  nHits_params_validator.add<string>("fit", "type");
786  nHits_params_validator.add<vector<double>>("nHits_params", {1., 1.});
787  std::vector<edm::ParameterSet> nHits_params_default;
788  edm::ParameterSet nHits_params1;
789  nHits_params1.addParameter<string>("fit", "hgcalEM");
790  nHits_params1.addParameter<vector<double>>("nHits_params", {157.522, 0.090});
791  nHits_params_default.push_back(nHits_params1);
792  edm::ParameterSet nHits_params2;
793  nHits_params2.addParameter<string>("fit", "hgcalHad");
794  nHits_params2.addParameter<vector<double>>("nHits_params", {159.295, 0.178});
795  nHits_params_default.push_back(nHits_params2);
796  edm::ParameterSet nHits_params3;
797  nHits_params3.addParameter<string>("fit", "hf");
798  nHits_params3.addParameter<vector<double>>("nHits_params", {165.706, 0.153});
799  nHits_params_default.push_back(nHits_params3);
800  desc.addVPSet("nHits_to_nvtx_params", nHits_params_validator, nHits_params_default);
801 
802  edm::ParameterSetDescription nvtx_params_validator;
803  nvtx_params_validator.add<string>("calo", "type");
804  nvtx_params_validator.add<string>("iEta", "etaregion");
805  nvtx_params_validator.add<vector<double>>("nvtx_params", {1., 1.});
806  std::vector<edm::ParameterSet> nvtx_params_default;
807  edm::ParameterSet nvtx_params1;
808  nvtx_params1.addParameter<string>("calo", "hgcalEM");
809  nvtx_params1.addParameter<string>("iEta", "er1p4to1p8");
810  nvtx_params1.addParameter<vector<double>>("nvtx_params", {-0.011772, 0.004142});
811  nvtx_params_default.push_back(nvtx_params1);
812  edm::ParameterSet nvtx_params2;
813  nvtx_params2.addParameter<string>("calo", "hgcalEM");
814  nvtx_params2.addParameter<string>("iEta", "er1p8to2p1");
815  nvtx_params2.addParameter<vector<double>>("nvtx_params", {-0.015488, 0.005410});
816  nvtx_params_default.push_back(nvtx_params2);
817  edm::ParameterSet nvtx_params3;
818  nvtx_params3.addParameter<string>("calo", "hgcalEM");
819  nvtx_params3.addParameter<string>("iEta", "er2p1to2p4");
820  nvtx_params3.addParameter<vector<double>>("nvtx_params", {-0.021150, 0.006078});
821  nvtx_params_default.push_back(nvtx_params3);
822  edm::ParameterSet nvtx_params4;
823  nvtx_params4.addParameter<string>("calo", "hgcalEM");
824  nvtx_params4.addParameter<string>("iEta", "er2p4to2p7");
825  nvtx_params4.addParameter<vector<double>>("nvtx_params", {-0.015705, 0.005339});
826  nvtx_params_default.push_back(nvtx_params4);
827  edm::ParameterSet nvtx_params5;
828  nvtx_params5.addParameter<string>("calo", "hgcalEM");
829  nvtx_params5.addParameter<string>("iEta", "er2p7to3p1");
830  nvtx_params5.addParameter<vector<double>>("nvtx_params", {-0.018492, 0.005620});
831  nvtx_params_default.push_back(nvtx_params5);
832  edm::ParameterSet nvtx_params6;
833  nvtx_params6.addParameter<string>("calo", "hgcalHad");
834  nvtx_params6.addParameter<string>("iEta", "er1p4to1p8");
835  nvtx_params6.addParameter<vector<double>>("nvtx_params", {0.005675, 0.000615});
836  nvtx_params_default.push_back(nvtx_params6);
837  edm::ParameterSet nvtx_params7;
838  nvtx_params7.addParameter<string>("calo", "hgcalHad");
839  nvtx_params7.addParameter<string>("iEta", "er1p8to2p1");
840  nvtx_params7.addParameter<vector<double>>("nvtx_params", {0.004560, 0.001099});
841  nvtx_params_default.push_back(nvtx_params7);
842  edm::ParameterSet nvtx_params8;
843  nvtx_params8.addParameter<string>("calo", "hgcalHad");
844  nvtx_params8.addParameter<string>("iEta", "er2p1to2p4");
845  nvtx_params8.addParameter<vector<double>>("nvtx_params", {0.000036, 0.001608});
846  nvtx_params_default.push_back(nvtx_params8);
847  edm::ParameterSet nvtx_params9;
848  nvtx_params9.addParameter<string>("calo", "hgcalHad");
849  nvtx_params9.addParameter<string>("iEta", "er2p4to2p7");
850  nvtx_params9.addParameter<vector<double>>("nvtx_params", {0.000869, 0.001754});
851  nvtx_params_default.push_back(nvtx_params9);
852  edm::ParameterSet nvtx_params10;
853  nvtx_params10.addParameter<string>("calo", "hgcalHad");
854  nvtx_params10.addParameter<string>("iEta", "er2p7to3p1");
855  nvtx_params10.addParameter<vector<double>>("nvtx_params", {-0.006574, 0.003134});
856  nvtx_params_default.push_back(nvtx_params10);
857  edm::ParameterSet nvtx_params11;
858  nvtx_params11.addParameter<string>("calo", "hf");
859  nvtx_params11.addParameter<string>("iEta", "er29to33");
860  nvtx_params11.addParameter<vector<double>>("nvtx_params", {-0.203291, 0.044096});
861  nvtx_params_default.push_back(nvtx_params11);
862  edm::ParameterSet nvtx_params12;
863  nvtx_params12.addParameter<string>("calo", "hf");
864  nvtx_params12.addParameter<string>("iEta", "er34to37");
865  nvtx_params12.addParameter<vector<double>>("nvtx_params", {-0.210922, 0.045628});
866  nvtx_params_default.push_back(nvtx_params12);
867  edm::ParameterSet nvtx_params13;
868  nvtx_params13.addParameter<string>("calo", "hf");
869  nvtx_params13.addParameter<string>("iEta", "er38to41");
870  nvtx_params13.addParameter<vector<double>>("nvtx_params", {-0.229562, 0.050560});
871  nvtx_params_default.push_back(nvtx_params13);
872  desc.addVPSet("nvtx_to_PU_sub_params", nvtx_params_validator, nvtx_params_default);
873 
874  desc.add<vector<double>>("jetPtBins", {0.0, 5.0, 7.5, 10.0, 12.5, 15.0, 17.5, 20.0, 22.5, 25.0, 27.5,
875  30.0, 35.0, 40.0, 45.0, 50.0, 55.0, 60.0, 65.0, 70.0, 75.0, 80.0,
876  85.0, 90.0, 95.0, 100.0, 110.0, 120.0, 130.0, 140.0, 150.0, 160.0, 170.0,
877  180.0, 190.0, 200.0, 225.0, 250.0, 275.0, 300.0, 325.0, 400.0, 500.0});
878  desc.add<vector<double>>("absEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
879  desc.add<vector<double>>(
880  "jetCalibrationsBarrel",
881  {2.459, 2.320, 2.239, 2.166, 2.100, 2.040, 1.986, 1.937, 1.892, 1.852, 1.816, 1.768, 1.714, 1.670, 1.633, 1.603,
882  1.578, 1.557, 1.540, 1.525, 1.513, 1.502, 1.493, 1.486, 1.479, 1.470, 1.460, 1.452, 1.445, 1.439, 1.433, 1.427,
883  1.422, 1.417, 1.411, 1.403, 1.390, 1.377, 1.365, 1.352, 1.327, 1.284, 4.695, 3.320, 2.751, 2.361, 2.093, 1.908,
884  1.781, 1.694, 1.633, 1.591, 1.562, 1.533, 1.511, 1.499, 1.492, 1.486, 1.482, 1.478, 1.474, 1.470, 1.467, 1.463,
885  1.459, 1.456, 1.452, 1.447, 1.440, 1.433, 1.425, 1.418, 1.411, 1.404, 1.397, 1.390, 1.382, 1.370, 1.352, 1.334,
886  1.316, 1.298, 1.262, 1.200, 5.100, 3.538, 2.892, 2.448, 2.143, 1.933, 1.789, 1.689, 1.620, 1.572, 1.539, 1.506,
887  1.482, 1.469, 1.460, 1.455, 1.450, 1.446, 1.442, 1.438, 1.434, 1.431, 1.427, 1.423, 1.420, 1.414, 1.407, 1.400,
888  1.392, 1.385, 1.378, 1.370, 1.363, 1.356, 1.348, 1.336, 1.317, 1.299, 1.281, 1.263, 1.226, 1.162, 3.850, 3.438,
889  3.211, 3.017, 2.851, 2.708, 2.585, 2.479, 2.388, 2.310, 2.243, 2.159, 2.072, 2.006, 1.956, 1.917, 1.887, 1.863,
890  1.844, 1.828, 1.814, 1.802, 1.791, 1.782, 1.773, 1.760, 1.744, 1.729, 1.714, 1.699, 1.685, 1.670, 1.656, 1.641,
891  1.627, 1.602, 1.566, 1.530, 1.494, 1.458, 1.386, 1.260});
892  desc.add<vector<double>>("absEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
893  desc.add<vector<double>>(
894  "jetCalibrationsHGCal",
895  {5.604, 4.578, 4.061, 3.647, 3.314, 3.047, 2.832, 2.660, 2.521, 2.410, 2.320, 2.216, 2.120, 2.056,
896  2.013, 1.983, 1.961, 1.945, 1.932, 1.922, 1.913, 1.905, 1.898, 1.891, 1.884, 1.874, 1.861, 1.848,
897  1.835, 1.822, 1.810, 1.797, 1.784, 1.771, 1.759, 1.736, 1.704, 1.673, 1.641, 1.609, 1.545, 1.434,
898  4.385, 3.584, 3.177, 2.849, 2.584, 2.370, 2.197, 2.057, 1.944, 1.853, 1.780, 1.695, 1.616, 1.564,
899  1.530, 1.507, 1.491, 1.480, 1.472, 1.466, 1.462, 1.459, 1.456, 1.453, 1.451, 1.447, 1.443, 1.439,
900  1.435, 1.431, 1.427, 1.423, 1.419, 1.416, 1.412, 1.405, 1.395, 1.385, 1.376, 1.366, 1.346, 1.312,
901  562.891, 68.647, 17.648, 5.241, 2.223, 1.490, 1.312, 1.270, 1.260, 1.259, 1.259, 1.260, 1.263, 1.265,
902  1.267, 1.269, 1.271, 1.273, 1.275, 1.277, 1.279, 1.281, 1.283, 1.285, 1.287, 1.290, 1.295, 1.299,
903  1.303, 1.307, 1.311, 1.315, 1.319, 1.323, 1.328, 1.335, 1.345, 1.355, 1.366, 1.376, 1.397, 1.433});
904  desc.add<vector<double>>("absEtaBinsHF", {3.00, 3.60, 6.00});
905  desc.add<vector<double>>(
906  "jetCalibrationsHF",
907  {8.169, 6.873, 6.155, 5.535, 5.001, 4.539, 4.141, 3.798, 3.501, 3.245, 3.024, 2.748, 2.463, 2.249,
908  2.090, 1.971, 1.881, 1.814, 1.763, 1.725, 1.695, 1.673, 1.655, 1.642, 1.631, 1.618, 1.605, 1.596,
909  1.588, 1.581, 1.575, 1.569, 1.563, 1.557, 1.551, 1.541, 1.527, 1.513, 1.498, 1.484, 1.456, 1.406,
910  2.788, 2.534, 2.388, 2.258, 2.141, 2.037, 1.945, 1.862, 1.788, 1.722, 1.664, 1.587, 1.503, 1.436,
911  1.382, 1.339, 1.305, 1.277, 1.255, 1.237, 1.223, 1.211, 1.201, 1.193, 1.186, 1.178, 1.170, 1.164,
912  1.159, 1.154, 1.151, 1.147, 1.144, 1.141, 1.138, 1.133, 1.126, 1.118, 1.111, 1.104, 1.090, 1.064});
913  desc.add<vector<double>>("tauPtBins", {0.0, 5.0, 7.5, 10.0, 12.5, 15.0, 20.0, 25.0, 30.0, 35.0,
914  40.0, 45.0, 50.0, 55.0, 60.0, 70.0, 80.0, 100.0, 150.0, 200.0});
915  desc.add<vector<double>>("tauAbsEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
916  desc.add<vector<double>>("tauCalibrationsBarrel",
917  {1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067,
918  1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106,
919  1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.102,
920  1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102,
921  1.102, 1.102, 1.102, 1.102, 1.102, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139});
922  desc.add<vector<double>>("tauAbsEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
923  desc.add<vector<double>>(
924  "tauCalibrationsHGCal",
925  {1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384,
926  1.384, 1.384, 1.384, 1.384, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473,
927  1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133,
928  1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133});
929 
930  descriptions.addWithDefaultLabel(desc);
931 }
932 
933 //define this as a plug-in
static float towerEta(int ieta)
Definition: CaloTools.cc:201
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::vector< double > jetPtBins
std::map< std::string, TF1 > hgcalEM_nvtx_to_PU_sub_funcs
std::map< std::string, TF1 > hgcalHad_nvtx_to_PU_sub_funcs
std::map< std::string, std::map< std::string, TF1 > > all_nvtx_to_PU_sub_funcs
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static constexpr int nBarrelPhi
void setTowerEt(float towerEtIn)
std::vector< double > tauPtBins
void setTowerIEta(int towerIEtaIn)
std::vector< double > jetCalibrationsBarrel
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::vector< std::vector< double > > calibrationsBarrel
static constexpr int nSTEta
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:70
std::map< std::string, TF1 > nHits_to_nvtx_funcs
static constexpr int nHgcalEta
delete x;
Definition: CaloConfig.h:22
std::vector< double > tauCalibrationsHGCal
static constexpr int nJets
std::vector< double > absEtaBinsHGCal
static const int kHFBegin
Definition: CaloTools.h:39
edm::EDGetTokenT< l1t::HGCalTowerBxCollection > hgcalTowerToken_
const_iterator begin(int bx) const
static constexpr int nHfPhi
std::vector< edm::ParameterSet > nHits_to_nvtx_params
BXVector< HGCalTower > HGCalTowerBxCollection
Definition: HGCalTower.h:10
float get_tau_pt_calibration(const float &tau_pt, const float &tau_eta) const
int iEvent
Definition: GenABIO.cc:224
void makeST_hf(const float hfTowers[nHfEta/2][nHfPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
std::vector< std::vector< double > > calibrationsHGCal
static constexpr int nSTPhi
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:136
static const int kHFEnd
Definition: CaloTools.h:40
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > jetCalibrationsHF
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< std::vector< double > > calibrationsHF
void setTauEt(float tauEtIn)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setTowerPhi(float towerPhiIn)
#define M_PI
void setJetIPhi(int jetIPhiIn)
unsigned int id
std::vector< double > tauAbsEtaBinsHGCal
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > decoderTag_
std::map< std::string, TF1 > hf_nvtx_to_PU_sub_funcs
std::vector< std::vector< double > > tauPtCalibrationsHGCal
void makeST(const float GCTintTowers[nBarrelEta/2][nBarrelPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
static constexpr int nHgcalPhi
void setTowerEta(float towerEtaIn)
const_iterator end(int bx) const
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void makeST_hgcal(const float hgcalTowers[nHgcalEta/2][nHgcalPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
edm::EDGetTokenT< HcalTrigPrimDigiCollection > hfToken_
std::vector< double > absEtaBinsBarrel
static constexpr int nBarrelEta
std::vector< double > tauCalibrationsBarrel
static constexpr int nHfEta
Phase2L1CaloJetEmulator(const edm::ParameterSet &)
std::vector< double > tauAbsEtaBinsBarrel
std::vector< edm::ParameterSet > nvtx_to_PU_sub_params
jetInfo getRegion(GCTsupertower_t temp[nSTEta][nSTPhi], float TTseedThreshold)
std::vector< double > absEtaBinsHF
void setTowerIPhi(int towerIPhiIn)
void setJetEt(float jetEtIn)
void setJetIEta(int jetIEtaIn)
std::vector< std::vector< double > > tauPtCalibrationsBarrel
edm::EDGetTokenT< l1tp2::CaloTowerCollection > caloTowerToken_
bool compareByEt(l1tp2::Phase2L1CaloJet i, l1tp2::Phase2L1CaloJet j)
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
void setJetPhi(float jetPhiIn)
Definition: Common.h:9
float get_jet_pt_calibration(const float &jet_pt, const float &jet_eta) const
std::vector< double > jetCalibrationsHGCal
void setJetEta(float jetEtaIn)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38