CMS 3D CMS Logo

L1TowerCalibrator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1CaloTrigger
4 // Class: L1TowerCalibrator
5 //
33 //
34 // Original Author: Tyler Ruggles
35 // Created: Thu Nov 15 2018
36 // $Id$
37 //
38 //
39 
47 
48 #include <iostream>
49 
54 
57 
59 
60 #include "TFile.h"
61 #include "TF1.h"
62 
64 public:
65  explicit L1TowerCalibrator(const edm::ParameterSet&);
66 
67 private:
68  void produce(edm::Event&, const edm::EventSetup&) override;
69 
70  const double HcalTpEtMin;
71  const double EcalTpEtMin;
72  const double HGCalHadTpEtMin;
73  const double HGCalEmTpEtMin;
74  const double HFTpEtMin;
75  const double puThreshold;
76  const double puThresholdL1eg;
77  const double puThresholdHcalMin;
78  const double puThresholdHcalMax;
79  const double puThresholdEcalMin;
80  const double puThresholdEcalMax;
81  const double puThresholdHGCalEMMin;
82  const double puThresholdHGCalEMMax;
83  const double puThresholdHGCalHadMin;
84  const double puThresholdHGCalHadMax;
85  const double puThresholdHFMin;
86  const double puThresholdHFMax;
87  const double barrelSF;
88  const double hgcalSF;
89  const double hfSF;
90  const bool debug;
91  const bool skipCalibrations;
92 
95 
99 
103 
104  // nHits to nvtx functions
105  std::vector<edm::ParameterSet> nHits_to_nvtx_params;
106  std::map<std::string, TF1> nHits_to_nvtx_funcs;
107 
108  // nvtx to PU subtraction functions
109  std::vector<edm::ParameterSet> nvtx_to_PU_sub_params;
110  std::map<std::string, TF1> ecal_nvtx_to_PU_sub_funcs;
111  std::map<std::string, TF1> hcal_nvtx_to_PU_sub_funcs;
112  std::map<std::string, TF1> hgcalEM_nvtx_to_PU_sub_funcs;
113  std::map<std::string, TF1> hgcalHad_nvtx_to_PU_sub_funcs;
114  std::map<std::string, TF1> hf_nvtx_to_PU_sub_funcs;
115  std::map<std::string, std::map<std::string, TF1> > all_nvtx_to_PU_sub_funcs;
116 };
117 
119  : HcalTpEtMin(iConfig.getParameter<double>("HcalTpEtMin")),
120  EcalTpEtMin(iConfig.getParameter<double>("EcalTpEtMin")),
121  HGCalHadTpEtMin(iConfig.getParameter<double>("HGCalHadTpEtMin")),
122  HGCalEmTpEtMin(iConfig.getParameter<double>("HGCalEmTpEtMin")),
123  HFTpEtMin(iConfig.getParameter<double>("HFTpEtMin")),
124  puThreshold(iConfig.getParameter<double>("puThreshold")),
125  puThresholdL1eg(iConfig.getParameter<double>("puThresholdL1eg")),
126  puThresholdHcalMin(iConfig.getParameter<double>("puThresholdHcalMin")),
127  puThresholdHcalMax(iConfig.getParameter<double>("puThresholdHcalMax")),
128  puThresholdEcalMin(iConfig.getParameter<double>("puThresholdEcalMin")),
129  puThresholdEcalMax(iConfig.getParameter<double>("puThresholdEcalMax")),
130  puThresholdHGCalEMMin(iConfig.getParameter<double>("puThresholdHGCalEMMin")),
131  puThresholdHGCalEMMax(iConfig.getParameter<double>("puThresholdHGCalEMMax")),
132  puThresholdHGCalHadMin(iConfig.getParameter<double>("puThresholdHGCalHadMin")),
133  puThresholdHGCalHadMax(iConfig.getParameter<double>("puThresholdHGCalHadMax")),
134  puThresholdHFMin(iConfig.getParameter<double>("puThresholdHFMin")),
135  puThresholdHFMax(iConfig.getParameter<double>("puThresholdHFMax")),
136  barrelSF(iConfig.getParameter<double>("barrelSF")),
137  hgcalSF(iConfig.getParameter<double>("hgcalSF")),
138  hfSF(iConfig.getParameter<double>("hfSF")),
139  debug(iConfig.getParameter<bool>("debug")),
140  skipCalibrations(iConfig.getParameter<bool>("skipCalibrations")),
141  l1TowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("l1CaloTowers"))),
142  hgcalTowersToken_(
143  consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("L1HgcalTowersInputTag"))),
144  hcalToken_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("hcalDigis"))),
145  decoderTag_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>(edm::ESInputTag("", ""))),
146  nHits_to_nvtx_params(iConfig.getParameter<std::vector<edm::ParameterSet> >("nHits_to_nvtx_params")),
147  nvtx_to_PU_sub_params(iConfig.getParameter<std::vector<edm::ParameterSet> >("nvtx_to_PU_sub_params")) {
148  // Initialize the nHits --> nvtx functions
149  for (uint i = 0; i < nHits_to_nvtx_params.size(); i++) {
151  std::string calo = pset->getParameter<std::string>("fit");
152  nHits_to_nvtx_funcs[calo.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
153  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(0, pset->getParameter<std::vector<double> >("params").at(0));
154  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(1, pset->getParameter<std::vector<double> >("params").at(1));
155 
156  if (debug) {
157  printf(
158  "nHits_to_nvtx_params[%i]\n \
159  fit: %s \n \
160  p1: %f \n \
161  p2 %f \n",
162  i,
163  calo.c_str(),
164  nHits_to_nvtx_funcs[calo.c_str()].GetParameter(0),
165  nHits_to_nvtx_funcs[calo.c_str()].GetParameter(1));
166  }
167  }
168 
169  // Initialize the nvtx --> PU subtraction functions
175 
176  for (uint i = 0; i < nvtx_to_PU_sub_params.size(); i++) {
178  std::string calo = pset->getParameter<std::string>("calo");
179  std::string iEta = pset->getParameter<std::string>("iEta");
180  double p1 = pset->getParameter<std::vector<double> >("params").at(0);
181  double p2 = pset->getParameter<std::vector<double> >("params").at(1);
182 
183  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
184  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(0, p1);
185  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(1, p2);
186 
187  if (debug) {
188  printf(
189  "nvtx_to_PU_sub_params[%i]\n \
190  sub detector: %s \n \
191  iEta: %s \n \
192  p1: %f \n \
193  p2 %f \n",
194  i,
195  calo.c_str(),
196  iEta.c_str(),
197  p1,
198  p2);
199  }
200  }
201 
202  // Our two outputs, calibrated towers and estimated nvtx for fun
203  produces<l1tp2::CaloTowerCollection>("L1CaloTowerCalibratedCollection");
204  produces<double>("EstimatedNvtx");
205 }
206 
208  // Estimated number of vertices used for calibration estimattion
209  std::unique_ptr<double> EstimatedNvtx(new double);
210  // Calibrated output collection
211  std::unique_ptr<l1tp2::CaloTowerCollection> L1CaloTowerCalibratedCollection(new l1tp2::CaloTowerCollection);
212 
213  // Load the ECAL+HCAL tower sums coming from L1EGammaCrystalsEmulatorProducer.cc
215 
216  // HGCal info
218  hgcalTowers = (*hgcalTowersHandle.product());
219 
220  // HF Tower info
221  iEvent.getByToken(hcalToken_, hcalTowerHandle);
222 
223  // Barrel ECAL (unclustered) and HCAL
224  for (auto& hit : *l1CaloTowerHandle.product()) {
225  l1tp2::CaloTower l1Hit;
226  l1Hit.setEcalTowerEt(hit.ecalTowerEt());
227  l1Hit.setHcalTowerEt(hit.hcalTowerEt());
228  l1Hit.setL1egTowerEt(hit.l1egTowerEt());
229  // Add min ET thresholds for tower ET
230  if (l1Hit.ecalTowerEt() < EcalTpEtMin)
231  l1Hit.setEcalTowerEt(0.0);
232  if (l1Hit.hcalTowerEt() < HcalTpEtMin)
233  l1Hit.setHcalTowerEt(0.0);
234  l1Hit.setTowerIEta(hit.towerIEta());
235  l1Hit.setTowerIPhi(hit.towerIPhi());
236  l1Hit.setTowerEta(hit.towerEta());
237  l1Hit.setTowerPhi(hit.towerPhi());
238  l1Hit.setIsBarrel(hit.isBarrel());
239  l1Hit.setNL1eg(hit.nL1eg());
240  l1Hit.setL1egTrkSS(hit.l1egTrkSS());
241  l1Hit.setL1egTrkIso(hit.l1egTrkIso());
242  l1Hit.setL1egStandaloneSS(hit.l1egStandaloneSS());
243  l1Hit.setL1egStandaloneIso(hit.l1egStandaloneIso());
244 
245  // FIXME There is an error in the L1EGammaCrystalsEmulatorProducer.cc which is
246  // returning towers with minimal ECAL energy, and no HCAL energy with these
247  // iEta/iPhi coordinates and eta = -88.653152 and phi = -99.000000.
248  // Skip these for the time being until the upstream code has been debugged
249  if ((int)l1Hit.towerIEta() == -1016 && (int)l1Hit.towerIPhi() == -962)
250  continue;
251 
252  (*L1CaloTowerCalibratedCollection).push_back(l1Hit);
253  if (debug)
254  printf("Barrel tower iEta %i iPhi %i eta %f phi %f ecal_et %f hcal_et_sum %f\n",
255  (int)l1Hit.towerIEta(),
256  (int)l1Hit.towerIPhi(),
257  l1Hit.towerEta(),
258  l1Hit.towerPhi(),
259  l1Hit.ecalTowerEt(),
260  l1Hit.hcalTowerEt());
261  }
262 
263  // Loop over HGCalTowers and create L1CaloTowers for them and add to collection
264  // This iterator is taken from the PF P2 group
265  // https://github.com/p2l1pfp/cmssw/blob/170808db68038d53794bc65fdc962f8fc337a24d/L1Trigger/Phase2L1ParticleFlow/plugins/L1TPFCaloProducer.cc#L278-L289
266  for (auto it = hgcalTowers.begin(0), ed = hgcalTowers.end(0); it != ed; ++it) {
267  // skip lowest ET towers
268  if (it->etEm() < HGCalEmTpEtMin && it->etHad() < HGCalHadTpEtMin)
269  continue;
270 
271  l1tp2::CaloTower l1Hit;
272  // Set energies normally, but need to zero if below threshold
273  if (it->etEm() < HGCalEmTpEtMin)
274  l1Hit.setEcalTowerEt(0.);
275  else
276  l1Hit.setEcalTowerEt(it->etEm());
277 
278  if (it->etHad() < HGCalHadTpEtMin)
279  l1Hit.setHcalTowerEt(0.);
280  else
281  l1Hit.setHcalTowerEt(it->etHad());
282 
283  l1Hit.setTowerEta(it->eta());
284  l1Hit.setTowerPhi(it->phi());
285  l1Hit.setTowerIEta(-98); // -98 mean HGCal
286  l1Hit.setTowerIPhi(-98);
287  l1Hit.setIsBarrel(false);
288  (*L1CaloTowerCalibratedCollection).push_back(l1Hit);
289  if (debug)
290  printf("HGCal tower iEta %i iPhi %i eta %f phi %f ecal_et %f hcal_et_sum %f\n",
291  (int)l1Hit.towerIEta(),
292  (int)l1Hit.towerIPhi(),
293  l1Hit.towerEta(),
294  l1Hit.towerPhi(),
295  l1Hit.ecalTowerEt(),
296  l1Hit.hcalTowerEt());
297  }
298 
299  // Loop over Hcal HF tower inputs and create L1CaloTowers and add to
300  // L1CaloTowerCalibratedCollection collection
301  const auto& decoder = iSetup.getData(decoderTag_);
302  for (const auto& hit : *hcalTowerHandle.product()) {
303  HcalTrigTowerDetId id = hit.id();
304  double et = decoder.hcaletValue(hit.id(), hit.t0());
305  if (et < HFTpEtMin)
306  continue;
307  // Only doing HF so skip outside range
308  if (abs(id.ieta()) < l1t::CaloTools::kHFBegin)
309  continue;
310  if (abs(id.ieta()) > l1t::CaloTools::kHFEnd)
311  continue;
312 
313  l1tp2::CaloTower l1Hit;
314  l1Hit.setEcalTowerEt(0.);
315  l1Hit.setHcalTowerEt(et);
317  l1Hit.setTowerPhi(l1t::CaloTools::towerPhi(id.ieta(), id.iphi()));
318  l1Hit.setTowerIEta(id.ieta());
319  l1Hit.setTowerIPhi(id.iphi());
320  l1Hit.setIsBarrel(false);
321  (*L1CaloTowerCalibratedCollection).push_back(l1Hit);
322 
323  if (debug)
324  printf("HCAL HF tower iEta %i iPhi %i eta %f phi %f ecal_et %f hcal_et_sum %f\n",
325  (int)l1Hit.towerIEta(),
326  (int)l1Hit.towerIPhi(),
327  l1Hit.towerEta(),
328  l1Hit.towerPhi(),
329  l1Hit.ecalTowerEt(),
330  l1Hit.hcalTowerEt());
331  }
332 
333  // N Tower totals
334  // For mapping to estimated nvtx in event
335  int i_ecal_hits_leq_threshold = 0;
336  int i_hgcalEM_hits_leq_threshold = 0;
337  int i_hcal_hits_leq_threshold = 0;
338  int i_hgcalHad_hits_leq_threshold = 0;
339  int i_hf_hits_leq_threshold = 0;
340 
341  // Loop over the collection containing all hits
342  // and calculate the number of hits falling into the
343  // "less than or equal" nTowers variable which maps to
344  // estimated number of vertices
345  for (auto& l1CaloTower : (*L1CaloTowerCalibratedCollection)) {
346  // Barrel ECAL
347  if (l1CaloTower.ecalTowerEt() > 0. && l1CaloTower.towerIEta() != -98) {
348  if (l1CaloTower.ecalTowerEt() <= puThresholdEcalMax && l1CaloTower.ecalTowerEt() >= puThresholdEcalMin) {
349  i_ecal_hits_leq_threshold++;
350  }
351  }
352 
353  // HGCal EM
354  if (l1CaloTower.ecalTowerEt() > 0. && l1CaloTower.towerIEta() == -98) {
355  if (l1CaloTower.ecalTowerEt() <= puThresholdHGCalEMMax && l1CaloTower.ecalTowerEt() >= puThresholdHGCalEMMin) {
356  i_hgcalEM_hits_leq_threshold++;
357  }
358  }
359 
360  // Barrel HCAL
361  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() != -98 &&
362  abs(l1CaloTower.towerEta()) < 2.0) // abs(eta) < 2 just keeps us out of HF
363  {
364  if (l1CaloTower.hcalTowerEt() <= puThresholdHcalMax && l1CaloTower.hcalTowerEt() >= puThresholdHcalMin) {
365  i_hcal_hits_leq_threshold++;
366  }
367  }
368 
369  // HGCal Had
370  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() == -98) {
371  if (l1CaloTower.hcalTowerEt() <= puThresholdHGCalHadMax && l1CaloTower.hcalTowerEt() >= puThresholdHGCalHadMin) {
372  i_hgcalHad_hits_leq_threshold++;
373  }
374  }
375 
376  // HF
377  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() != -98 &&
378  abs(l1CaloTower.towerEta()) > 2.0) // abs(eta) > 2 keeps us out of barrel HF
379  {
380  if (l1CaloTower.hcalTowerEt() <= puThresholdHFMax && l1CaloTower.hcalTowerEt() >= puThresholdHFMin) {
381  i_hf_hits_leq_threshold++;
382  }
383  }
384  }
385 
386  // For each subdetector, map to the estimated number of vertices
387  double ecal_nvtx = nHits_to_nvtx_funcs["ecal"].Eval(i_ecal_hits_leq_threshold);
388  double hcal_nvtx = nHits_to_nvtx_funcs["hcal"].Eval(i_hcal_hits_leq_threshold);
389  double hgcalEM_nvtx = nHits_to_nvtx_funcs["hgcalEM"].Eval(i_hgcalEM_hits_leq_threshold);
390  double hgcalHad_nvtx = nHits_to_nvtx_funcs["hgcalHad"].Eval(i_hgcalHad_hits_leq_threshold);
391  double hf_nvtx = nHits_to_nvtx_funcs["hf"].Eval(i_hf_hits_leq_threshold);
392  // Make sure all values are >= 0
393  if (ecal_nvtx < 0)
394  ecal_nvtx = 0;
395  if (hcal_nvtx < 0)
396  hcal_nvtx = 0;
397  if (hgcalEM_nvtx < 0)
398  hgcalEM_nvtx = 0;
399  if (hgcalHad_nvtx < 0)
400  hgcalHad_nvtx = 0;
401  if (hf_nvtx < 0)
402  hf_nvtx = 0;
403  // Best estimate is avg of all except HF.
404  // This is b/c HF currently has such poor prediction power, it only degrades the avg result
405  // NEW, with 10_3_X, hgcal and HF has the best results based on the values I took...
406  // skip ECAL and HCAL
407  //*EstimatedNvtx = ( ecal_nvtx + hcal_nvtx + hgcalEM_nvtx + hgcalHad_nvtx + hf_nvtx ) / 3.;
408  *EstimatedNvtx = (hgcalEM_nvtx + hgcalHad_nvtx + hf_nvtx) / 3.;
409 
410  if (debug) {
411  double lumi = iEvent.eventAuxiliary().luminosityBlock();
412  double event = iEvent.eventAuxiliary().event();
413 
414  printf(
415  "L1TowerCalibrater: lumi %.0f evt %.0f nTowers for subdetecters \
416  \nECAL: %i --> nvtx = %.1f \
417  \nHGCal EM: %i --> nvtx = %.1f \
418  \nHCAL: %i --> nvtx = %.1f \
419  \nHGCal Had: %i --> nvtx = %.1f \
420  \nHCAL HF: %i --> nvtx = %.1f \
421  \nEstimated Nvtx = %.1f\n",
422  lumi,
423  event,
424  i_ecal_hits_leq_threshold,
425  ecal_nvtx,
426  i_hgcalEM_hits_leq_threshold,
427  hgcalEM_nvtx,
428  i_hcal_hits_leq_threshold,
429  hcal_nvtx,
430  i_hgcalHad_hits_leq_threshold,
431  hgcalHad_nvtx,
432  i_hf_hits_leq_threshold,
433  hf_nvtx,
434  *EstimatedNvtx);
435  }
436 
437  // Use estimated number of vertices to subtract off PU contributions
438  // to each and every hit. In cases where the energy would go negative,
439  // limit this to zero.
440  if (!skipCalibrations) // skipCalibrations simply passes the towers through
441  {
442  for (auto& l1CaloTower : (*L1CaloTowerCalibratedCollection)) {
443  // Barrel ECAL eta slices
444  if (l1CaloTower.ecalTowerEt() > 0. && l1CaloTower.towerIEta() != -98) {
445  if (abs(l1CaloTower.towerIEta()) <= 3) {
446  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
447  all_nvtx_to_PU_sub_funcs["ecal"]["er1to3"].Eval(*EstimatedNvtx) * barrelSF);
448  }
449  if (abs(l1CaloTower.towerIEta()) <= 6 && abs(l1CaloTower.towerIEta()) >= 4) {
450  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
451  all_nvtx_to_PU_sub_funcs["ecal"]["er4to6"].Eval(*EstimatedNvtx) * barrelSF);
452  }
453  if (abs(l1CaloTower.towerIEta()) <= 9 && abs(l1CaloTower.towerIEta()) >= 7) {
454  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
455  all_nvtx_to_PU_sub_funcs["ecal"]["er7to9"].Eval(*EstimatedNvtx) * barrelSF);
456  }
457  if (abs(l1CaloTower.towerIEta()) <= 12 && abs(l1CaloTower.towerIEta()) >= 10) {
458  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
459  all_nvtx_to_PU_sub_funcs["ecal"]["er10to12"].Eval(*EstimatedNvtx) * barrelSF);
460  }
461  if (abs(l1CaloTower.towerIEta()) <= 15 && abs(l1CaloTower.towerIEta()) >= 13) {
462  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
463  all_nvtx_to_PU_sub_funcs["ecal"]["er13to15"].Eval(*EstimatedNvtx) * barrelSF);
464  }
465  if (abs(l1CaloTower.towerIEta()) <= 18 && abs(l1CaloTower.towerIEta()) >= 16) {
466  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
467  all_nvtx_to_PU_sub_funcs["ecal"]["er16to18"].Eval(*EstimatedNvtx) * barrelSF);
468  }
469  }
470 
471  // HGCal EM eta slices
472  if (l1CaloTower.ecalTowerEt() > 0. && l1CaloTower.towerIEta() == -98) {
473  if (abs(l1CaloTower.towerEta()) <= 1.8) {
474  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
475  all_nvtx_to_PU_sub_funcs["hgcalEM"]["er1p4to1p8"].Eval(*EstimatedNvtx) * hgcalSF);
476  }
477  if (abs(l1CaloTower.towerEta()) <= 2.1 && abs(l1CaloTower.towerEta()) > 1.8) {
478  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
479  all_nvtx_to_PU_sub_funcs["hgcalEM"]["er1p8to2p1"].Eval(*EstimatedNvtx) * hgcalSF);
480  }
481  if (abs(l1CaloTower.towerEta()) <= 2.4 && abs(l1CaloTower.towerEta()) > 2.1) {
482  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
483  all_nvtx_to_PU_sub_funcs["hgcalEM"]["er2p1to2p4"].Eval(*EstimatedNvtx) * hgcalSF);
484  }
485  if (abs(l1CaloTower.towerEta()) <= 2.7 && abs(l1CaloTower.towerEta()) > 2.4) {
486  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
487  all_nvtx_to_PU_sub_funcs["hgcalEM"]["er2p4to2p7"].Eval(*EstimatedNvtx) * hgcalSF);
488  }
489  if (abs(l1CaloTower.towerEta()) <= 3.1 && abs(l1CaloTower.towerEta()) > 2.7) {
490  l1CaloTower.setEcalTowerEt(l1CaloTower.ecalTowerEt() -
491  all_nvtx_to_PU_sub_funcs["hgcalEM"]["er2p7to3p1"].Eval(*EstimatedNvtx) * hgcalSF);
492  }
493  }
494 
495  // Barrel HCAL eta slices
496  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() != -98 &&
497  abs(l1CaloTower.towerEta()) < 2.0) // abs(eta) < 2 just keeps us out of HF
498  {
499  if (abs(l1CaloTower.towerIEta()) <= 3) {
500  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
501  all_nvtx_to_PU_sub_funcs["hcal"]["er1to3"].Eval(*EstimatedNvtx) * barrelSF);
502  }
503  if (abs(l1CaloTower.towerIEta()) <= 6 && abs(l1CaloTower.towerIEta()) >= 4) {
504  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
505  all_nvtx_to_PU_sub_funcs["hcal"]["er4to6"].Eval(*EstimatedNvtx) * barrelSF);
506  }
507  if (abs(l1CaloTower.towerIEta()) <= 9 && abs(l1CaloTower.towerIEta()) >= 7) {
508  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
509  all_nvtx_to_PU_sub_funcs["hcal"]["er7to9"].Eval(*EstimatedNvtx) * barrelSF);
510  }
511  if (abs(l1CaloTower.towerIEta()) <= 12 && abs(l1CaloTower.towerIEta()) >= 10) {
512  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
513  all_nvtx_to_PU_sub_funcs["hcal"]["er10to12"].Eval(*EstimatedNvtx) * barrelSF);
514  }
515  if (abs(l1CaloTower.towerIEta()) <= 15 && abs(l1CaloTower.towerIEta()) >= 13) {
516  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
517  all_nvtx_to_PU_sub_funcs["hcal"]["er13to15"].Eval(*EstimatedNvtx) * barrelSF);
518  }
519  if (abs(l1CaloTower.towerIEta()) <= 18 && abs(l1CaloTower.towerIEta()) >= 16) {
520  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
521  all_nvtx_to_PU_sub_funcs["hcal"]["er16to18"].Eval(*EstimatedNvtx) * barrelSF);
522  }
523  }
524 
525  // HGCal Had eta slices
526  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() == -98) {
527  if (abs(l1CaloTower.towerEta()) <= 1.8) {
528  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
529  all_nvtx_to_PU_sub_funcs["hgcalHad"]["er1p4to1p8"].Eval(*EstimatedNvtx) * hgcalSF);
530  }
531  if (abs(l1CaloTower.towerEta()) <= 2.1 && abs(l1CaloTower.towerEta()) > 1.8) {
532  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
533  all_nvtx_to_PU_sub_funcs["hgcalHad"]["er1p8to2p1"].Eval(*EstimatedNvtx) * hgcalSF);
534  }
535  if (abs(l1CaloTower.towerEta()) <= 2.4 && abs(l1CaloTower.towerEta()) > 2.1) {
536  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
537  all_nvtx_to_PU_sub_funcs["hgcalHad"]["er2p1to2p4"].Eval(*EstimatedNvtx) * hgcalSF);
538  }
539  if (abs(l1CaloTower.towerEta()) <= 2.7 && abs(l1CaloTower.towerEta()) > 2.4) {
540  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
541  all_nvtx_to_PU_sub_funcs["hgcalHad"]["er2p4to2p7"].Eval(*EstimatedNvtx) * hgcalSF);
542  }
543  if (abs(l1CaloTower.towerEta()) <= 3.1 && abs(l1CaloTower.towerEta()) > 2.7) {
544  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
545  all_nvtx_to_PU_sub_funcs["hgcalHad"]["er2p7to3p1"].Eval(*EstimatedNvtx) * hgcalSF);
546  }
547  }
548 
549  // HF eta slices
550  if (l1CaloTower.hcalTowerEt() > 0. && l1CaloTower.towerIEta() != -98 &&
551  abs(l1CaloTower.towerEta()) > 2.0) // abs(eta) > 2 keeps us out of barrel HF
552  {
553  if (abs(l1CaloTower.towerIEta()) <= 33 && abs(l1CaloTower.towerIEta()) >= 29) {
554  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
555  all_nvtx_to_PU_sub_funcs["hf"]["er29to33"].Eval(*EstimatedNvtx) * hfSF);
556  }
557  if (abs(l1CaloTower.towerIEta()) <= 37 && abs(l1CaloTower.towerIEta()) >= 34) {
558  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
559  all_nvtx_to_PU_sub_funcs["hf"]["er34to37"].Eval(*EstimatedNvtx) * hfSF);
560  }
561  if (abs(l1CaloTower.towerIEta()) <= 41 && abs(l1CaloTower.towerIEta()) >= 38) {
562  l1CaloTower.setHcalTowerEt(l1CaloTower.hcalTowerEt() -
563  all_nvtx_to_PU_sub_funcs["hf"]["er38to41"].Eval(*EstimatedNvtx) * hfSF);
564  }
565  }
566 
567  // Make sure none are negative
568  if (l1CaloTower.ecalTowerEt() < 0.)
569  l1CaloTower.setEcalTowerEt(0.);
570  if (l1CaloTower.hcalTowerEt() < 0.)
571  l1CaloTower.setHcalTowerEt(0.);
572  }
573  }
574 
575  iEvent.put(std::move(EstimatedNvtx), "EstimatedNvtx");
576  iEvent.put(std::move(L1CaloTowerCalibratedCollection), "L1CaloTowerCalibratedCollection");
577 }
578 
static float towerEta(int ieta)
Definition: CaloTools.cc:201
const double puThresholdHGCalEMMin
const double HGCalEmTpEtMin
float ecalTowerEt() const
Definition: CaloTower.h:35
const double puThresholdHcalMax
const double HFTpEtMin
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const double puThresholdEcalMax
static float towerPhi(int ieta, int iphi)
Definition: CaloTools.cc:208
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const double HGCalHadTpEtMin
std::map< std::string, TF1 > nHits_to_nvtx_funcs
std::map< std::string, TF1 > hf_nvtx_to_PU_sub_funcs
const double puThresholdL1eg
const bool skipCalibrations
void setEcalTowerEt(float et)
Definition: CaloTower.h:49
int towerIPhi() const
Definition: CaloTower.h:37
void setTowerIPhi(int iPhi)
Definition: CaloTower.h:51
T const * product() const
Definition: Handle.h:70
std::vector< CaloTower > CaloTowerCollection
Definition: CaloTower.h:83
void setL1egStandaloneIso(int staIso)
Definition: CaloTower.h:60
const edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > decoderTag_
delete x;
Definition: CaloConfig.h:22
const double puThresholdEcalMin
void setL1egTrkIso(int trkIso)
Definition: CaloTower.h:58
const edm::EDGetTokenT< l1t::HGCalTowerBxCollection > hgcalTowersToken_
const double puThresholdHFMin
static const int kHFBegin
Definition: CaloTools.h:39
const double HcalTpEtMin
const_iterator begin(int bx) const
std::map< std::string, TF1 > hgcalEM_nvtx_to_PU_sub_funcs
BXVector< HGCalTower > HGCalTowerBxCollection
Definition: HGCalTower.h:10
std::vector< edm::ParameterSet > nHits_to_nvtx_params
void setHcalTowerEt(float et)
Definition: CaloTower.h:50
void setTowerIEta(int iEta)
Definition: CaloTower.h:52
L1TowerCalibrator(const edm::ParameterSet &)
const double EcalTpEtMin
int iEvent
Definition: GenABIO.cc:224
void setTowerPhi(float phi)
Definition: CaloTower.h:53
std::map< std::string, std::map< std::string, TF1 > > all_nvtx_to_PU_sub_funcs
const double puThresholdHGCalEMMax
static const int kHFEnd
Definition: CaloTools.h:40
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< std::string, TF1 > hcal_nvtx_to_PU_sub_funcs
int towerIEta() const
Definition: CaloTower.h:38
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setL1egTrkSS(int trkSS)
Definition: CaloTower.h:57
void setIsBarrel(bool isBarrel)
Definition: CaloTower.h:61
void setTowerEta(float eta)
Definition: CaloTower.h:54
unsigned int id
const edm::EDGetTokenT< l1tp2::CaloTowerCollection > l1TowerToken_
const double puThresholdHcalMin
std::vector< edm::ParameterSet > nvtx_to_PU_sub_params
#define debug
Definition: HDRShower.cc:19
float towerPhi() const
Definition: CaloTower.h:39
void produce(edm::Event &, const edm::EventSetup &) override
edm::Handle< l1tp2::CaloTowerCollection > l1CaloTowerHandle
float hcalTowerEt() const
Definition: CaloTower.h:36
void setL1egTowerEt(float et)
Definition: CaloTower.h:55
edm::Handle< l1t::HGCalTowerBxCollection > hgcalTowersHandle
std::map< std::string, TF1 > hgcalHad_nvtx_to_PU_sub_funcs
const_iterator end(int bx) const
HLT enums.
const double puThresholdHFMax
const double puThresholdHGCalHadMax
const edm::EDGetTokenT< HcalTrigPrimDigiCollection > hcalToken_
const double puThreshold
const double puThresholdHGCalHadMin
edm::Handle< HcalTrigPrimDigiCollection > hcalTowerHandle
std::map< std::string, TF1 > ecal_nvtx_to_PU_sub_funcs
void setNL1eg(int n)
Definition: CaloTower.h:56
float towerEta() const
Definition: CaloTower.h:40
l1t::HGCalTowerBxCollection hgcalTowers
void setL1egStandaloneSS(int staSS)
Definition: CaloTower.h:59
def move(src, dest)
Definition: eostools.py:511
Definition: Common.h:9
Definition: event.py:1