CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
JetIDHelper.cc
Go to the documentation of this file.
2 
4 
9 
12 // JetIDHelper needs a much more detailed description that the one in HcalTopology,
13 // so to be consistent, all needed constants are hardwired in JetIDHelper.cc itself
16 
17 #include "TMath.h"
18 #include <vector>
19 #include <numeric>
20 #include <iostream>
21 
22 using namespace std;
23 
24 // defining stuff useful only for this class (not needed in header)
25 namespace reco {
26 
27  namespace helper {
28 
29  // select2nd exists only in some std and boost implementations, so let's control our own fate
30  // and it can't be a non-static member function.
31  static double select2nd(std::map<int, double>::value_type const &pair) { return pair.second; }
32 
34 
36  return i.E > j.E;
37  }
38  std::atomic<int> JetIDHelper::sanity_checks_left_{100};
39  } // namespace helper
40 } // namespace reco
41 
43  useRecHits_ = pset.getParameter<bool>("useRecHits");
44  if (useRecHits_) {
45  hbheRecHitsColl_ = pset.getParameter<edm::InputTag>("hbheRecHitsColl");
46  hoRecHitsColl_ = pset.getParameter<edm::InputTag>("hoRecHitsColl");
47  hfRecHitsColl_ = pset.getParameter<edm::InputTag>("hfRecHitsColl");
48  ebRecHitsColl_ = pset.getParameter<edm::InputTag>("ebRecHitsColl");
49  eeRecHitsColl_ = pset.getParameter<edm::InputTag>("eeRecHitsColl");
50  }
51  initValues();
52 
53  input_HBHERecHits_token_ = iC.consumes<HBHERecHitCollection>(hbheRecHitsColl_);
54  input_HORecHits_token_ = iC.consumes<HORecHitCollection>(hoRecHitsColl_);
55  input_HFRecHits_token_ = iC.consumes<HFRecHitCollection>(hfRecHitsColl_);
56  input_EBRecHits_token_ = iC.consumes<EBRecHitCollection>(ebRecHitsColl_);
57  input_EERecHits_token_ = iC.consumes<EERecHitCollection>(eeRecHitsColl_);
58 
59  hcal_topo_token_ = iC.esConsumes();
60 }
61 
63  fHPD_ = -1.0;
64  fRBX_ = -1.0;
65  n90Hits_ = -1;
66  fSubDetector1_ = -1.0;
67  fSubDetector2_ = -1.0;
68  fSubDetector3_ = -1.0;
69  fSubDetector4_ = -1.0;
70  restrictedEMF_ = -1.0;
71  nHCALTowers_ = -1;
72  nECALTowers_ = -1;
73  approximatefHPD_ = -1.0;
74  approximatefRBX_ = -1.0;
75  hitsInN90_ = -1;
76  fEB_ = fEE_ = fHB_ = fHE_ = fHO_ = fLong_ = fShort_ = -1.0;
77  fLS_ = fHFOOT_ = -1.0;
78 }
79 
81  iDesc
82  .ifValue(edm::ParameterDescription<bool>("useRecHits", true, true),
83  true >> (edm::ParameterDescription<edm::InputTag>("hbheRecHitsColl", edm::InputTag(), true) and
84  edm::ParameterDescription<edm::InputTag>("hoRecHitsColl", edm::InputTag(), true) and
85  edm::ParameterDescription<edm::InputTag>("hfRecHitsColl", edm::InputTag(), true) and
86  edm::ParameterDescription<edm::InputTag>("ebRecHitsColl", edm::InputTag(), true) and
88  ->setComment(
89  "If using RecHits to calculate the precise jet ID variables that need them, "
90  "their sources need to be specified");
91 }
92 
94  const edm::EventSetup &setup,
95  const reco::CaloJet &jet,
96  const int iDbg) {
97  initValues();
98 
99  // --------------------------------------------------
100  // 1) jet ID variable derived from existing fractions
101  // --------------------------------------------------
102 
103  double E_EM = TMath::Max(float(0.), jet.emEnergyInHF()) + TMath::Max(float(0.), jet.emEnergyInEB()) +
104  TMath::Max(float(0.), jet.emEnergyInEE());
105  double E_Had = TMath::Max(float(0.), jet.hadEnergyInHB()) + TMath::Max(float(0.), jet.hadEnergyInHE()) +
106  TMath::Max(float(0.), jet.hadEnergyInHO()) + TMath::Max(float(0.), jet.hadEnergyInHF());
107  if (E_Had + E_EM > 0)
108  restrictedEMF_ = E_EM / (E_EM + E_Had);
109  if (iDbg > 1)
110  cout << "jet pT: " << jet.pt() << ", eT: " << jet.et() << ", E: " << jet.energy() << " rEMF: " << restrictedEMF_
111  << endl;
112 
113  // ------------------------
114  // 2) tower based variables
115  // ------------------------
116  vector<subtower> subtowers, Ecal_subtowers, Hcal_subtowers, HO_subtowers;
117  vector<double> HPD_energies, RBX_energies;
118 
119  classifyJetTowers(
120  event, jet, subtowers, Ecal_subtowers, Hcal_subtowers, HO_subtowers, HPD_energies, RBX_energies, iDbg);
121  if (iDbg > 1) {
122  cout << "E:";
123  for (unsigned int i = 0; i < subtowers.size(); ++i)
124  cout << " " << subtowers[i].E << "," << subtowers[i].Nhit;
125  cout << "\nECal_E:";
126  for (unsigned int i = 0; i < Ecal_subtowers.size(); ++i)
127  cout << " " << Ecal_subtowers[i].E << "," << Ecal_subtowers[i].Nhit;
128  cout << "\nHCal_E:";
129  for (unsigned int i = 0; i < Hcal_subtowers.size(); ++i)
130  cout << " " << Hcal_subtowers[i].E << "," << Hcal_subtowers[i].Nhit;
131  cout << "\nHO_E:";
132  for (unsigned int i = 0; i < HO_subtowers.size(); ++i)
133  cout << " " << HO_subtowers[i].E << "," << HO_subtowers[i].Nhit;
134  cout << "\nHPD_E:";
135  for (unsigned int i = 0; i < HPD_energies.size(); ++i)
136  cout << " " << HPD_energies[i];
137  cout << "\nRBX_E:";
138  for (unsigned int i = 0; i < RBX_energies.size(); ++i)
139  cout << " " << RBX_energies[i];
140  cout << endl;
141  }
142 
143  // counts
144  hitsInN90_ = hitsInNCarrying(0.9, subtowers);
145  nHCALTowers_ = Hcal_subtowers.size();
146  vector<subtower>::const_iterator it;
147  it = find_if(Ecal_subtowers.begin(), Ecal_subtowers.end(), hasNonPositiveE);
148  nECALTowers_ = it - Ecal_subtowers.begin(); // ignores negative energies from HF!
149 
150  // energy fractions
151  double max_HPD_energy = 0., max_RBX_energy = 0.;
152  std::vector<double>::const_iterator it_max_HPD_energy = std::max_element(HPD_energies.begin(), HPD_energies.end());
153  std::vector<double>::const_iterator it_max_RBX_energy = std::max_element(RBX_energies.begin(), RBX_energies.end());
154  if (it_max_HPD_energy != HPD_energies.end()) {
155  max_HPD_energy = *it_max_HPD_energy;
156  }
157  if (it_max_RBX_energy != RBX_energies.end()) {
158  max_RBX_energy = *it_max_RBX_energy;
159  }
160  if (jet.energy() > 0) {
161  if (!HPD_energies.empty())
162  approximatefHPD_ = max_HPD_energy / jet.energy();
163  if (!RBX_energies.empty())
164  approximatefRBX_ = max_HPD_energy / jet.energy();
165  }
166 
167  // -----------------------
168  // 3) cell based variables
169  // -----------------------
170  if (useRecHits_) {
171  vector<double> energies, subdet_energies, Ecal_energies, Hcal_energies, HO_energies;
172  double LS_bad_energy, HF_OOT_energy;
173  classifyJetComponents(event,
174  setup,
175  jet,
176  energies,
177  subdet_energies,
178  Ecal_energies,
179  Hcal_energies,
180  HO_energies,
181  HPD_energies,
182  RBX_energies,
183  LS_bad_energy,
184  HF_OOT_energy,
185  iDbg);
186 
187  // counts
188  n90Hits_ = nCarrying(0.9, energies);
189 
190  // energy fractions
191  if (jet.energy() > 0) {
192  if (!HPD_energies.empty())
193  fHPD_ = max_HPD_energy / jet.energy();
194  if (!RBX_energies.empty())
195  fRBX_ = max_RBX_energy / jet.energy();
196  if (!subdet_energies.empty())
197  fSubDetector1_ = subdet_energies.at(0) / jet.energy();
198  if (subdet_energies.size() > 1)
199  fSubDetector2_ = subdet_energies.at(1) / jet.energy();
200  if (subdet_energies.size() > 2)
201  fSubDetector3_ = subdet_energies.at(2) / jet.energy();
202  if (subdet_energies.size() > 3)
203  fSubDetector4_ = subdet_energies.at(3) / jet.energy();
204  fLS_ = LS_bad_energy / jet.energy();
205  fHFOOT_ = HF_OOT_energy / jet.energy();
206 
207  if (iDbg > 0 && sanity_checks_left_.load(std::memory_order_acquire) > 0) {
208  --sanity_checks_left_;
209  double EH_sum = accumulate(Ecal_energies.begin(), Ecal_energies.end(), 0.);
210  EH_sum = accumulate(Hcal_energies.begin(), Hcal_energies.end(), EH_sum);
211  double EHO_sum = accumulate(HO_energies.begin(), HO_energies.end(), EH_sum);
212  if (jet.energy() > 0.001 && abs(EH_sum / jet.energy() - 1) > 0.01 && abs(EHO_sum / jet.energy() - 1) > 0.01)
213  edm::LogWarning("BadInput") << "Jet energy (" << jet.energy()
214  << ") does not match the total energy in the recHits (" << EHO_sum
215  << ", or without HO: " << EH_sum << ") . Are these the right recHits? "
216  << "Were jet energy corrections mistakenly applied before jet ID? A bug?";
217  if (iDbg > 1)
218  cout << "Sanity check - E: " << jet.energy() << " =? EH_sum: " << EH_sum << " / EHO_sum: " << EHO_sum << endl;
219  }
220  }
221 
222  if (iDbg > 1) {
223  cout << "DBG - fHPD: " << fHPD_ << ", fRBX: " << fRBX_ << ", nh90: " << n90Hits_ << ", fLS: " << fLS_
224  << ", fHFOOT: " << fHFOOT_ << endl;
225  cout << " -~fHPD: " << approximatefHPD_ << ", ~fRBX: " << approximatefRBX_ << ", hits in n90: " << hitsInN90_
226  << endl;
227  cout << " - nHCALTowers: " << nHCALTowers_ << ", nECALTowers: " << nECALTowers_
228  << "; subdet fractions: " << fSubDetector1_ << ", " << fSubDetector2_ << ", " << fSubDetector3_ << ", "
229  << fSubDetector4_ << endl;
230  }
231  }
232 }
233 
234 unsigned int reco::helper::JetIDHelper::nCarrying(double fraction, const std::vector<double> &descending_energies) {
235  double totalE = 0;
236  for (unsigned int i = 0; i < descending_energies.size(); ++i)
237  totalE += descending_energies[i];
238 
239  double runningE = 0;
240  unsigned int NC = descending_energies.size();
241 
242  // slightly odd loop structure avoids round-off problems when runningE never catches up with totalE
243  for (unsigned int i = descending_energies.size(); i > 0; --i) {
244  runningE += descending_energies[i - 1];
245  if (runningE < (1 - fraction) * totalE)
246  NC = i - 1;
247  }
248  return NC;
249 }
250 
252  const std::vector<subtower> &descending_towers) {
253  double totalE = 0;
254  for (unsigned int i = 0; i < descending_towers.size(); ++i)
255  totalE += descending_towers[i].E;
256 
257  double runningE = 0;
258  unsigned int NH = 0;
259 
260  // slightly odd loop structure avoids round-off problems when runningE never catches up with totalE
261  for (unsigned int i = descending_towers.size(); i > 0; --i) {
262  runningE += descending_towers[i - 1].E;
263  if (runningE >= (1 - fraction) * totalE)
264  NH += descending_towers[i - 1].Nhit;
265  }
266  return NH;
267 }
268 
270  const edm::EventSetup &setup,
271  const reco::CaloJet &jet,
272  vector<double> &energies,
273  vector<double> &subdet_energies,
274  vector<double> &Ecal_energies,
275  vector<double> &Hcal_energies,
276  vector<double> &HO_energies,
277  vector<double> &HPD_energies,
278  vector<double> &RBX_energies,
279  double &LS_bad_energy,
280  double &HF_OOT_energy,
281  const int iDbg) {
282  energies.clear();
283  subdet_energies.clear();
284  Ecal_energies.clear();
285  Hcal_energies.clear();
286  HO_energies.clear();
287  HPD_energies.clear();
288  RBX_energies.clear();
289  LS_bad_energy = HF_OOT_energy = 0.;
290 
291  edm::ESHandle<HcalTopology> theHcalTopology = setup.getHandle(hcal_topo_token_);
292 
293  std::map<int, double> HPD_energy_map, RBX_energy_map;
294  vector<double> EB_energies, EE_energies, HB_energies, HE_energies, short_energies, long_energies;
300  // the jet only contains DetIds, so first read recHit collection
301  event.getByToken(input_HBHERecHits_token_, HBHERecHits);
302  event.getByToken(input_HORecHits_token_, HORecHits);
303  event.getByToken(input_HFRecHits_token_, HFRecHits);
304  event.getByToken(input_EBRecHits_token_, EBRecHits);
305  event.getByToken(input_EERecHits_token_, EERecHits);
306  if (iDbg > 2)
307  cout << "# of rechits found - HBHE: " << HBHERecHits->size() << ", HO: " << HORecHits->size()
308  << ", HF: " << HFRecHits->size() << ", EB: " << EBRecHits->size() << ", EE: " << EERecHits->size() << endl;
309 
310  vector<CaloTowerPtr> towers = jet.getCaloConstituents();
311  int nTowers = towers.size();
312  if (iDbg > 9)
313  cout << "In classifyJetComponents. # of towers found: " << nTowers << endl;
314 
315  for (int iTower = 0; iTower < nTowers; iTower++) {
316  CaloTowerPtr &tower = towers[iTower];
317 
318  int nCells = tower->constituentsSize();
319  if (iDbg)
320  cout << "tower #" << iTower << " has " << nCells << " cells. "
321  << "It's at iEta: " << tower->ieta() << ", iPhi: " << tower->iphi() << endl;
322 
323  const vector<DetId> &cellIDs = tower->constituents(); // cell == recHit
324 
325  for (int iCell = 0; iCell < nCells; ++iCell) {
326  DetId::Detector detNum = cellIDs[iCell].det();
327  if (detNum == DetId::Hcal) {
328  HcalDetId HcalID = cellIDs[iCell];
329  HcalSubdetector HcalNum = HcalID.subdet();
330  double hitE = 0;
331  if (HcalNum == HcalOuter) {
332  HORecHitCollection::const_iterator theRecHit = HORecHits->find(HcalID);
333  if (theRecHit == HORecHits->end()) {
334  edm::LogWarning("UnexpectedEventContents") << "Can't find the HO recHit with ID: " << HcalID;
335  continue;
336  }
337  hitE = theRecHit->energy();
338  HO_energies.push_back(hitE);
339 
340  } else if (HcalNum == HcalForward) {
341  HFRecHitCollection::const_iterator theRecHit = HFRecHits->find(HcalID);
342  if (theRecHit == HFRecHits->end()) {
343  edm::LogWarning("UnexpectedEventContents") << "Can't find the HF recHit with ID: " << HcalID;
344  continue;
345  }
346  hitE = theRecHit->energy();
347  if (iDbg > 4)
348  cout << "hit #" << iCell << " is HF , E: " << hitE << " iEta: " << theRecHit->id().ieta()
349  << ", depth: " << theRecHit->id().depth() << ", iPhi: " << theRecHit->id().iphi();
350 
351  if (HcalID.depth() == 1)
352  long_energies.push_back(hitE);
353  else
354  short_energies.push_back(hitE);
355 
356  uint32_t flags = theRecHit->flags();
357  if (flags & (1 << HcalCaloFlagLabels::HFLongShort))
358  LS_bad_energy += hitE;
362  HF_OOT_energy += hitE;
363  if (iDbg > 4 && flags)
364  cout << "flags: " << flags << " -> LS_bad_energy: " << LS_bad_energy << ", HF_OOT_energy: " << HF_OOT_energy
365  << endl;
366 
367  } else { // HBHE
368 
369  HBHERecHitCollection::const_iterator theRecHit = HBHERecHits->find(HcalID);
370  if (theRecHit == HBHERecHits->end()) {
371  edm::LogWarning("UnexpectedEventContents") << "Can't find the HBHE recHit with ID: " << HcalID;
372  continue;
373  }
374  hitE = theRecHit->energy();
375  int iEta = theRecHit->id().ieta();
376  int depth = theRecHit->id().depth();
377  Region region = HBHE_region(theRecHit->id().rawId());
378  int hitIPhi = theRecHit->id().iphi();
379  if (iDbg > 3)
380  cout << "hit #" << iCell << " is HBHE, E: " << hitE << " iEta: " << iEta << ", depth: " << depth
381  << ", iPhi: " << theRecHit->id().iphi() << " -> " << region;
382 
383  if (theHcalTopology->mergedDepth29(theRecHit->id()))
384  hitE /= 2; // Depth 3 at the HE forward edge is split over tower 28 & 29, and jet reco. assigns half each
385 
386  int iHPD = 100 * region;
387  int iRBX = 100 * region + ((hitIPhi + 1) % 72) / 4; // 71,72,1,2 are in the same RBX module
388 
389  if (std::abs(iEta) >= 21) {
390  if ((0x1 & hitIPhi) == 0) {
391  edm::LogError("CodeAssumptionsViolated") << "Bug?! Jet ID code assumes no even iPhi recHits at HE edges";
392  return;
393  }
394  bool oddnessIEta = HBHE_oddness(iEta, depth);
395  bool upperIPhi = ((hitIPhi % 4) == 1 || (hitIPhi % 4) == 2); // the upper iPhi indices in the HE wedge
396  // remap the iPhi so it fits the one in the inner HE regions, change in needed in two cases:
397  // 1) in the upper iPhis of the module, the even iEtas belong to the higher iPhi
398  // 2) in the loewr iPhis of the module, the odd iEtas belong to the higher iPhi
399  if (upperIPhi != oddnessIEta)
400  ++hitIPhi;
401  // note that hitIPhi could not be 72 before, so it's still in the legal range [1,72]
402  }
403  iHPD += hitIPhi;
404 
405  // book the energies
406  HPD_energy_map[iHPD] += hitE;
407  RBX_energy_map[iRBX] += hitE;
408  if (iDbg > 5)
409  cout << " --> H[" << iHPD << "]=" << HPD_energy_map[iHPD] << ", R[" << iRBX << "]=" << RBX_energy_map[iRBX];
410  if (iDbg > 1)
411  cout << endl;
412 
413  if (region == HBneg || region == HBpos)
414  HB_energies.push_back(hitE);
415  else
416  HE_energies.push_back(hitE);
417 
418  } // if HBHE
419  if (hitE == 0)
420  edm::LogWarning("UnexpectedEventContents") << "HCal hitE==0? (or unknown subdetector?)";
421 
422  } // if HCAL
423 
424  else if (detNum == DetId::Ecal) {
425  int EcalNum = cellIDs[iCell].subdetId();
426  double hitE = 0;
427  if (EcalNum == 1) {
428  EBDetId EcalID = cellIDs[iCell];
429  EBRecHitCollection::const_iterator theRecHit = EBRecHits->find(EcalID);
430  if (theRecHit == EBRecHits->end()) {
431  edm::LogWarning("UnexpectedEventContents") << "Can't find the EB recHit with ID: " << EcalID;
432  continue;
433  }
434  hitE = theRecHit->energy();
435  EB_energies.push_back(hitE);
436  } else if (EcalNum == 2) {
437  EEDetId EcalID = cellIDs[iCell];
438  EERecHitCollection::const_iterator theRecHit = EERecHits->find(EcalID);
439  if (theRecHit == EERecHits->end()) {
440  edm::LogWarning("UnexpectedEventContents") << "Can't find the EE recHit with ID: " << EcalID;
441  continue;
442  }
443  hitE = theRecHit->energy();
444  EE_energies.push_back(hitE);
445  }
446  if (hitE == 0)
447  edm::LogWarning("UnexpectedEventContents") << "ECal hitE==0? (or unknown subdetector?)";
448  if (iDbg > 6)
449  cout << "EcalNum: " << EcalNum << " hitE: " << hitE << endl;
450  } //
451  } // loop on cells
452  } // loop on towers
453 
454  /* Disabling check until HO is accounted for in EMF. Check was used in CMSSW_2, where HE was excluded.
455  double expHcalE = jet.energy() * (1-jet.emEnergyFraction());
456  if( totHcalE + expHcalE > 0 &&
457  TMath::Abs( totHcalE - expHcalE ) > 0.01 &&
458  ( totHcalE - expHcalE ) / ( totHcalE + expHcalE ) > 0.0001 ) {
459  edm::LogWarning("CodeAssumptionsViolated")<<"failed to account for all Hcal energies"
460  <<totHcalE<<"!="<<expHcalE;
461  } */
462  // concatenate Hcal and Ecal lists
463  Hcal_energies.insert(Hcal_energies.end(), HB_energies.begin(), HB_energies.end());
464  Hcal_energies.insert(Hcal_energies.end(), HE_energies.begin(), HE_energies.end());
465  Hcal_energies.insert(Hcal_energies.end(), short_energies.begin(), short_energies.end());
466  Hcal_energies.insert(Hcal_energies.end(), long_energies.begin(), long_energies.end());
467  Ecal_energies.insert(Ecal_energies.end(), EB_energies.begin(), EB_energies.end());
468  Ecal_energies.insert(Ecal_energies.end(), EE_energies.begin(), EE_energies.end());
469 
470  // sort the energies
471  // std::sort( Hcal_energies.begin(), Hcal_energies.end(), greater<double>() );
472  // std::sort( Ecal_energies.begin(), Ecal_energies.end(), greater<double>() );
473 
474  // put the energy sums (the 2nd entry in each pair of the maps) into the output vectors and sort them
476  HPD_energy_map.begin(), HPD_energy_map.end(), std::inserter(HPD_energies, HPD_energies.end()), select2nd);
477  // std::select2nd<std::map<int,double>::value_type>());
478  // std::sort( HPD_energies.begin(), HPD_energies.end(), greater<double>() );
480  RBX_energy_map.begin(), RBX_energy_map.end(), std::inserter(RBX_energies, RBX_energies.end()), select2nd);
481  // std::select2nd<std::map<int,double>::value_type>());
482  // std::sort( RBX_energies.begin(), RBX_energies.end(), greater<double>() );
483 
484  energies.insert(energies.end(), Hcal_energies.begin(), Hcal_energies.end());
485  energies.insert(energies.end(), Ecal_energies.begin(), Ecal_energies.end());
486  energies.insert(energies.end(), HO_energies.begin(), HO_energies.end());
487  std::sort(energies.begin(), energies.end(), greater<double>());
488 
489  // prepare sub detector energies, then turn them into fractions
490  fEB_ = std::accumulate(EB_energies.begin(), EB_energies.end(), 0.);
491  fEE_ = std::accumulate(EE_energies.begin(), EE_energies.end(), 0.);
492  fHB_ = std::accumulate(HB_energies.begin(), HB_energies.end(), 0.);
493  fHE_ = std::accumulate(HE_energies.begin(), HE_energies.end(), 0.);
494  fHO_ = std::accumulate(HO_energies.begin(), HO_energies.end(), 0.);
495  fShort_ = std::accumulate(short_energies.begin(), short_energies.end(), 0.);
496  fLong_ = std::accumulate(long_energies.begin(), long_energies.end(), 0.);
497  subdet_energies.push_back(fEB_);
498  subdet_energies.push_back(fEE_);
499  subdet_energies.push_back(fHB_);
500  subdet_energies.push_back(fHE_);
501  subdet_energies.push_back(fHO_);
502  subdet_energies.push_back(fShort_);
503  subdet_energies.push_back(fLong_);
504  std::sort(subdet_energies.begin(), subdet_energies.end(), greater<double>());
505  if (jet.energy() > 0) {
506  fEB_ /= jet.energy();
507  fEE_ /= jet.energy();
508  fHB_ /= jet.energy();
509  fHE_ /= jet.energy();
510  fHO_ /= jet.energy();
511  fShort_ /= jet.energy();
512  fLong_ /= jet.energy();
513  } else {
514  if (fEB_ > 0 || fEE_ > 0 || fHB_ > 0 || fHE_ > 0 || fHO_ > 0 || fShort_ > 0 || fLong_ > 0)
515  edm::LogError("UnexpectedEventContents") << "Jet ID Helper found energy in subdetectors and jet E <= 0";
516  }
517 }
518 
520  const reco::CaloJet &jet,
521  vector<subtower> &subtowers,
522  vector<subtower> &Ecal_subtowers,
523  vector<subtower> &Hcal_subtowers,
524  vector<subtower> &HO_subtowers,
525  vector<double> &HPD_energies,
526  vector<double> &RBX_energies,
527  const int iDbg) {
528  subtowers.clear();
529  Ecal_subtowers.clear();
530  Hcal_subtowers.clear();
531  HO_subtowers.clear();
532  HPD_energies.clear();
533  RBX_energies.clear();
534 
535  std::map<int, double> HPD_energy_map, RBX_energy_map;
536 
537  vector<CaloTowerPtr> towers = jet.getCaloConstituents();
538  int nTowers = towers.size();
539  if (iDbg > 9)
540  cout << "classifyJetTowers started. # of towers found: " << nTowers << endl;
541 
542  for (int iTower = 0; iTower < nTowers; iTower++) {
543  CaloTowerPtr &tower = towers[iTower];
544 
545  int nEM = 0, nHad = 0, nHO = 0;
546  const vector<DetId> &cellIDs = tower->constituents(); // cell == recHit
547  int nCells = cellIDs.size();
548  if (iDbg)
549  cout << "tower #" << iTower << " has " << nCells << " cells. "
550  << "It's at iEta: " << tower->ieta() << ", iPhi: " << tower->iphi() << endl;
551 
552  for (int iCell = 0; iCell < nCells; ++iCell) {
553  DetId::Detector detNum = cellIDs[iCell].det();
554  if (detNum == DetId::Hcal) {
555  HcalDetId HcalID = cellIDs[iCell];
556  HcalSubdetector HcalNum = HcalID.subdet();
557  if (HcalNum == HcalOuter) {
558  ++nHO;
559  } else {
560  ++nHad;
561  }
562  } else if (detNum == DetId::Ecal) {
563  ++nEM;
564  }
565  }
566 
567  double E_em = tower->emEnergy();
568  if (E_em != 0)
569  Ecal_subtowers.push_back(subtower(E_em, nEM));
570 
571  double E_HO = tower->outerEnergy();
572  if (E_HO != 0)
573  HO_subtowers.push_back(subtower(E_HO, nHO));
574 
575  double E_had = tower->hadEnergy();
576  if (E_had != 0) {
577  Hcal_subtowers.push_back(subtower(E_had, nHad));
578  // totHcalE += E_had;
579 
580  int iEta = tower->ieta();
581  Region reg = region(iEta);
582  int iPhi = tower->iphi();
583  if (iDbg > 3)
584  cout << "tower has E_had: " << E_had << " iEta: " << iEta << ", iPhi: " << iPhi << " -> " << reg;
585 
586  if (reg == HEneg || reg == HBneg || reg == HBpos || reg == HEpos) {
587  int oddnessIEta = HBHE_oddness(iEta);
588  if (oddnessIEta < 0)
589  break; // can't assign this tower to a single readout component
590 
591  int iHPD = 100 * reg;
592  int iRBX = 100 * reg + ((iPhi + 1) % 72) / 4; // 71,72,1,2 are in the same RBX module
593 
594  if ((reg == HEneg || reg == HEpos) && std::abs(iEta) >= 21) { // at low-granularity edge of HE
595  if ((0x1 & iPhi) == 0) {
596  edm::LogError("CodeAssumptionsViolated") << "Bug?! Jet ID code assumes no even iPhi recHits at HE edges";
597  return;
598  }
599  bool boolOddnessIEta = oddnessIEta;
600  bool upperIPhi = ((iPhi % 4) == 1 || (iPhi % 4) == 2); // the upper iPhi indices in the HE wedge
601  // remap the iPhi so it fits the one in the inner HE regions, change in needed in two cases:
602  // 1) in the upper iPhis of the module, the even IEtas belong to the higher iPhi
603  // 2) in the loewr iPhis of the module, the odd IEtas belong to the higher iPhi
604  if (upperIPhi != boolOddnessIEta)
605  ++iPhi;
606  // note that iPhi could not be 72 before, so it's still in the legal range [1,72]
607  } // if at low-granularity edge of HE
608  iHPD += iPhi;
609 
610  // book the energies
611  HPD_energy_map[iHPD] += E_had;
612  RBX_energy_map[iRBX] += E_had;
613  if (iDbg > 5)
614  cout << " --> H[" << iHPD << "]=" << HPD_energy_map[iHPD] << ", R[" << iRBX << "]=" << RBX_energy_map[iRBX];
615  } // HBHE
616  } // E_had > 0
617  } // loop on towers
618 
619  // sort the subtowers
620  // std::sort( Hcal_subtowers.begin(), Hcal_subtowers.end(), subtower_has_greater_E );
621  // std::sort( Ecal_subtowers.begin(), Ecal_subtowers.end(), subtower_has_greater_E );
622 
623  // put the energy sums (the 2nd entry in each pair of the maps) into the output vectors and sort them
625  HPD_energy_map.begin(), HPD_energy_map.end(), std::inserter(HPD_energies, HPD_energies.end()), select2nd);
626  // std::select2nd<std::map<int,double>::value_type>());
627  // std::sort( HPD_energies.begin(), HPD_energies.end(), greater<double>() );
629  RBX_energy_map.begin(), RBX_energy_map.end(), std::inserter(RBX_energies, RBX_energies.end()), select2nd);
630  // std::select2nd<std::map<int,double>::value_type>());
631  // std::sort( RBX_energies.begin(), RBX_energies.end(), greater<double>() );
632 
633  subtowers.insert(subtowers.end(), Hcal_subtowers.begin(), Hcal_subtowers.end());
634  subtowers.insert(subtowers.end(), Ecal_subtowers.begin(), Ecal_subtowers.end());
635  subtowers.insert(subtowers.end(), HO_subtowers.begin(), HO_subtowers.end());
636  std::sort(subtowers.begin(), subtowers.end(), subtower_has_greater_E);
637 }
638 
639 // ---------------------------------------------------------------------------------------------------------
640 // private helper functions to figure out detector & readout geometry as needed
641 
643  int ae = TMath::Abs(iEta);
644  if (ae == 29 && depth == 1)
645  ae += 1; // observed that: hits are at depths 1 & 2; 1 goes with the even pattern
646  return ae & 0x1;
647 }
648 
650  HcalDetId id(rawid);
651  if (id.subdet() == HcalEndcap) {
652  if (id.ieta() < 0)
653  return HEneg;
654  else
655  return HEpos;
656  }
657  if (id.subdet() == HcalBarrel && id.ieta() < 0)
658  return HBneg;
659  return HBpos;
660 }
661 
663  int ae = TMath::Abs(iEta);
664  if (ae == 29)
665  return -1; // can't figure it out without RecHits
666  return ae & 0x1;
667 }
668 
670  if (iEta == 16 || iEta == -16)
671  return unknown_region; // both HB and HE cells belong to these towers
672  if (iEta == 29 || iEta == -29)
673  return unknown_region; // both HE and HF cells belong to these towers
674  if (iEta <= -30)
675  return HFneg;
676  if (iEta >= 30)
677  return HFpos;
678  if (iEta <= -17)
679  return HEneg;
680  if (iEta >= 17)
681  return HEpos;
682  if (iEta < 0)
683  return HBneg;
684  return HBpos;
685 }
float hadEnergyInHE() const
Definition: CaloJet.h:103
float emEnergyInEE() const
Definition: CaloJet.h:109
void setComment(std::string const &value)
ParameterDescriptionNode * ifValue(ParameterDescription< T > const &switchParameter, std::unique_ptr< ParameterDescriptionCases< T >> cases)
Region region(int iEta)
Definition: JetIDHelper.cc:669
bool subtower_has_greater_E(reco::helper::JetIDHelper::subtower i, reco::helper::JetIDHelper::subtower j)
Definition: JetIDHelper.cc:35
void classifyJetComponents(const edm::Event &event, const edm::EventSetup &setup, const reco::CaloJet &jet, std::vector< double > &energies, std::vector< double > &subdet_energies, std::vector< double > &Ecal_energies, std::vector< double > &Hcal_energies, std::vector< double > &HO_energies, std::vector< double > &HPD_energies, std::vector< double > &RBX_energies, double &LS_bad_energy, double &HF_OOT_energy, const int iDbg=0)
Definition: JetIDHelper.cc:269
Jets made from CaloTowers.
Definition: CaloJet.h:27
double pt() const final
transverse momentum
uint16_t *__restrict__ id
std::vector< T >::const_iterator const_iterator
float emEnergyInHF() const
Definition: CaloJet.h:111
unsigned int nCarrying(double fraction, const std::vector< double > &descending_energies)
Definition: JetIDHelper.cc:234
Log< level::Error, false > LogError
virtual std::vector< CaloTowerPtr > getCaloConstituents() const
get all constituents
Definition: CaloJet.cc:78
void calculate(const edm::Event &event, const edm::EventSetup &setup, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:93
bool hasNonPositiveE(reco::helper::JetIDHelper::subtower x)
Definition: JetIDHelper.cc:33
Region HBHE_region(uint32_t)
Definition: JetIDHelper.cc:649
float hadEnergyInHO() const
Definition: CaloJet.h:101
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
float emEnergyInEB() const
Definition: CaloJet.h:107
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static double select2nd(std::map< int, double >::value_type const &pair)
Definition: JetIDHelper.cc:31
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
void fillDescription(edm::ParameterSetDescription &iDesc)
Definition: JetIDHelper.cc:80
void classifyJetTowers(const edm::Event &event, const reco::CaloJet &jet, std::vector< subtower > &subtowers, std::vector< subtower > &Ecal_subtowers, std::vector< subtower > &Hcal_subtowers, std::vector< subtower > &HO_subtowers, std::vector< double > &HPD_energies, std::vector< double > &RBX_energies, const int iDbg=0)
Definition: JetIDHelper.cc:519
Detector
Definition: DetId.h:24
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int HBHE_oddness(int iEta, int depth)
Definition: JetIDHelper.cc:642
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ nCells
double et() const final
transverse energy
float hadEnergyInHB() const
Definition: CaloJet.h:99
tuple cout
Definition: gather_cfg.py:144
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
Log< level::Warning, false > LogWarning
float hadEnergyInHF() const
Definition: CaloJet.h:105
unsigned int hitsInNCarrying(double fraction, const std::vector< subtower > &descending_towers)
Definition: JetIDHelper.cc:251
double energy() const final
energy
unsigned transform(const HcalDetId &id, unsigned transformCode)