CMS 3D CMS Logo

CaloTowersCreationAlgo.cc
Go to the documentation of this file.
2 
10 #include "Math/Interpolator.h"
11 
12 #include <cmath>
13 
14 //#define EDM_ML_DEBUG
15 
17  : theEBthreshold(-1000.),
18  theEEthreshold(-1000.),
19 
20  theUseEtEBTresholdFlag(false),
21  theUseEtEETresholdFlag(false),
22  theUseSymEBTresholdFlag(false),
23  theUseSymEETresholdFlag(false),
24 
25  theHcalThreshold(-1000.),
26  theHBthreshold(-1000.),
27  theHBthreshold1(-1000.),
28  theHBthreshold2(-1000.),
29  theHESthreshold(-1000.),
30  theHESthreshold1(-1000.),
31  theHEDthreshold(-1000.),
32  theHEDthreshold1(-1000.),
33  theHOthreshold0(-1000.),
34  theHOthresholdPlus1(-1000.),
35  theHOthresholdMinus1(-1000.),
36  theHOthresholdPlus2(-1000.),
37  theHOthresholdMinus2(-1000.),
38  theHF1threshold(-1000.),
39  theHF2threshold(-1000.),
40  theEBGrid(std::vector<double>(5, 10.)),
41  theEBWeights(std::vector<double>(5, 1.)),
42  theEEGrid(std::vector<double>(5, 10.)),
43  theEEWeights(std::vector<double>(5, 1.)),
44  theHBGrid(std::vector<double>(5, 10.)),
45  theHBWeights(std::vector<double>(5, 1.)),
46  theHESGrid(std::vector<double>(5, 10.)),
47  theHESWeights(std::vector<double>(5, 1.)),
48  theHEDGrid(std::vector<double>(5, 10.)),
49  theHEDWeights(std::vector<double>(5, 1.)),
50  theHOGrid(std::vector<double>(5, 10.)),
51  theHOWeights(std::vector<double>(5, 1.)),
52  theHF1Grid(std::vector<double>(5, 10.)),
53  theHF1Weights(std::vector<double>(5, 1.)),
54  theHF2Grid(std::vector<double>(5, 10.)),
55  theHF2Weights(std::vector<double>(5, 1.)),
56  theEBweight(1.),
57  theEEweight(1.),
58  theHBweight(1.),
59  theHESweight(1.),
60  theHEDweight(1.),
61  theHOweight(1.),
62  theHF1weight(1.),
63  theHF2weight(1.),
64  theEcutTower(-1000.),
65  theEBSumThreshold(-1000.),
66  theEESumThreshold(-1000.),
67  theEBEScale(50.),
68  theEEEScale(50.),
69  theHBEScale(50.),
70  theHESEScale(50.),
71  theHEDEScale(50.),
72  theHOEScale(50.),
73  theHF1EScale(50.),
74  theHF2EScale(50.),
75  theHcalTopology(nullptr),
76  theGeometry(nullptr),
77  theTowerConstituentsMap(nullptr),
78  theHcalAcceptSeverityLevel(0),
79  theRecoveredHcalHitsAreUsed(false),
80  theRecoveredEcalHitsAreUsed(false),
81  useRejectedHitsOnly(false),
82  theHcalAcceptSeverityLevelForRejectedHit(0),
83  useRejectedRecoveredHcalHits(0),
84  useRejectedRecoveredEcalHits(0),
86  theHOIsUsed(true),
87  // (for momentum reconstruction algorithm)
88  theMomConstrMethod(0),
89  theMomHBDepth(0.),
90  theMomHEDepth(0.),
91  theMomEBDepth(0.),
92  theMomEEDepth(0.),
93  theHcalPhase(0) {}
94 
96  double EEthreshold,
97 
98  bool useEtEBTreshold,
99  bool useEtEETreshold,
100  bool useSymEBTreshold,
101  bool useSymEETreshold,
102 
103  double HcalThreshold,
104  double HBthreshold,
105  double HBthreshold1,
106  double HBthreshold2,
107  double HESthreshold,
108  double HESthreshold1,
109  double HEDthreshold,
110  double HEDthreshold1,
111  double HOthreshold0,
112  double HOthresholdPlus1,
113  double HOthresholdMinus1,
114  double HOthresholdPlus2,
115  double HOthresholdMinus2,
116  double HF1threshold,
117  double HF2threshold,
118  double EBweight,
119  double EEweight,
120  double HBweight,
121  double HESweight,
122  double HEDweight,
123  double HOweight,
124  double HF1weight,
125  double HF2weight,
126  double EcutTower,
127  double EBSumThreshold,
128  double EESumThreshold,
129  bool useHO,
130  // (momentum reconstruction algorithm)
131  int momConstrMethod,
132  double momHBDepth,
133  double momHEDepth,
134  double momEBDepth,
135  double momEEDepth,
136  int hcalPhase)
137 
138  : theEBthreshold(EBthreshold),
139  theEEthreshold(EEthreshold),
140 
141  theUseEtEBTresholdFlag(useEtEBTreshold),
142  theUseEtEETresholdFlag(useEtEETreshold),
143  theUseSymEBTresholdFlag(useSymEBTreshold),
144  theUseSymEETresholdFlag(useSymEETreshold),
145 
146  theHcalThreshold(HcalThreshold),
147  theHBthreshold(HBthreshold),
148  theHBthreshold1(HBthreshold1),
149  theHBthreshold2(HBthreshold2),
150  theHESthreshold(HESthreshold),
151  theHESthreshold1(HESthreshold1),
152  theHEDthreshold(HEDthreshold),
153  theHEDthreshold1(HEDthreshold1),
154  theHOthreshold0(HOthreshold0),
155  theHOthresholdPlus1(HOthresholdPlus1),
156  theHOthresholdMinus1(HOthresholdMinus1),
157  theHOthresholdPlus2(HOthresholdPlus2),
158  theHOthresholdMinus2(HOthresholdMinus2),
159  theHF1threshold(HF1threshold),
160  theHF2threshold(HF2threshold),
161  theEBGrid(std::vector<double>(5, 10.)),
162  theEBWeights(std::vector<double>(5, 1.)),
163  theEEGrid(std::vector<double>(5, 10.)),
164  theEEWeights(std::vector<double>(5, 1.)),
165  theHBGrid(std::vector<double>(5, 10.)),
166  theHBWeights(std::vector<double>(5, 1.)),
167  theHESGrid(std::vector<double>(5, 10.)),
168  theHESWeights(std::vector<double>(5, 1.)),
169  theHEDGrid(std::vector<double>(5, 10.)),
170  theHEDWeights(std::vector<double>(5, 1.)),
171  theHOGrid(std::vector<double>(5, 10.)),
172  theHOWeights(std::vector<double>(5, 1.)),
173  theHF1Grid(std::vector<double>(5, 10.)),
174  theHF1Weights(std::vector<double>(5, 1.)),
175  theHF2Grid(std::vector<double>(5, 10.)),
176  theHF2Weights(std::vector<double>(5, 1.)),
177  theEBweight(EBweight),
178  theEEweight(EEweight),
179  theHBweight(HBweight),
180  theHESweight(HESweight),
181  theHEDweight(HEDweight),
182  theHOweight(HOweight),
183  theHF1weight(HF1weight),
184  theHF2weight(HF2weight),
185  theEcutTower(EcutTower),
186  theEBSumThreshold(EBSumThreshold),
187  theEESumThreshold(EESumThreshold),
188  theEBEScale(50.),
189  theEEEScale(50.),
190  theHBEScale(50.),
191  theHESEScale(50.),
192  theHEDEScale(50.),
193  theHOEScale(50.),
194  theHF1EScale(50.),
195  theHF2EScale(50.),
196  theHcalTopology(nullptr),
197  theGeometry(nullptr),
198  theTowerConstituentsMap(nullptr),
199  theHcalAcceptSeverityLevel(0),
200  theRecoveredHcalHitsAreUsed(false),
201  theRecoveredEcalHitsAreUsed(false),
202  useRejectedHitsOnly(false),
203  theHcalAcceptSeverityLevelForRejectedHit(0),
204  useRejectedRecoveredHcalHits(0),
205  useRejectedRecoveredEcalHits(0),
207  theHOIsUsed(useHO),
208  // (momentum reconstruction algorithm)
209  theMomConstrMethod(momConstrMethod),
210  theMomHBDepth(momHBDepth),
211  theMomHEDepth(momHEDepth),
212  theMomEBDepth(momEBDepth),
213  theMomEEDepth(momEEDepth),
214  theHcalPhase(hcalPhase) {}
215 
217  double EEthreshold,
218  bool useEtEBTreshold,
219  bool useEtEETreshold,
220  bool useSymEBTreshold,
221  bool useSymEETreshold,
222 
223  double HcalThreshold,
224  double HBthreshold,
225  double HBthreshold1,
226  double HBthreshold2,
227  double HESthreshold,
228  double HESthreshold1,
229  double HEDthreshold,
230  double HEDthreshold1,
231  double HOthreshold0,
232  double HOthresholdPlus1,
233  double HOthresholdMinus1,
234  double HOthresholdPlus2,
235  double HOthresholdMinus2,
236  double HF1threshold,
237  double HF2threshold,
238  const std::vector<double>& EBGrid,
239  const std::vector<double>& EBWeights,
240  const std::vector<double>& EEGrid,
241  const std::vector<double>& EEWeights,
242  const std::vector<double>& HBGrid,
243  const std::vector<double>& HBWeights,
244  const std::vector<double>& HESGrid,
245  const std::vector<double>& HESWeights,
246  const std::vector<double>& HEDGrid,
247  const std::vector<double>& HEDWeights,
248  const std::vector<double>& HOGrid,
249  const std::vector<double>& HOWeights,
250  const std::vector<double>& HF1Grid,
251  const std::vector<double>& HF1Weights,
252  const std::vector<double>& HF2Grid,
253  const std::vector<double>& HF2Weights,
254  double EBweight,
255  double EEweight,
256  double HBweight,
257  double HESweight,
258  double HEDweight,
259  double HOweight,
260  double HF1weight,
261  double HF2weight,
262  double EcutTower,
263  double EBSumThreshold,
264  double EESumThreshold,
265  bool useHO,
266  // (for the momentum construction algorithm)
267  int momConstrMethod,
268  double momHBDepth,
269  double momHEDepth,
270  double momEBDepth,
271  double momEEDepth,
272  int hcalPhase)
273 
274  : theEBthreshold(EBthreshold),
275  theEEthreshold(EEthreshold),
276 
277  theUseEtEBTresholdFlag(useEtEBTreshold),
278  theUseEtEETresholdFlag(useEtEETreshold),
279  theUseSymEBTresholdFlag(useSymEBTreshold),
280  theUseSymEETresholdFlag(useSymEETreshold),
281 
282  theHcalThreshold(HcalThreshold),
283  theHBthreshold(HBthreshold),
284  theHBthreshold1(HBthreshold1),
285  theHBthreshold2(HBthreshold2),
286  theHESthreshold(HESthreshold),
287  theHESthreshold1(HESthreshold1),
288  theHEDthreshold(HEDthreshold),
289  theHEDthreshold1(HEDthreshold1),
290  theHOthreshold0(HOthreshold0),
291  theHOthresholdPlus1(HOthresholdPlus1),
292  theHOthresholdMinus1(HOthresholdMinus1),
293  theHOthresholdPlus2(HOthresholdPlus2),
294  theHOthresholdMinus2(HOthresholdMinus2),
295  theHF1threshold(HF1threshold),
296  theHF2threshold(HF2threshold),
297  theEBGrid(EBGrid),
298  theEBWeights(EBWeights),
299  theEEGrid(EEGrid),
300  theEEWeights(EEWeights),
301  theHBGrid(HBGrid),
302  theHBWeights(HBWeights),
303  theHESGrid(HESGrid),
304  theHESWeights(HESWeights),
305  theHEDGrid(HEDGrid),
306  theHEDWeights(HEDWeights),
307  theHOGrid(HOGrid),
308  theHOWeights(HOWeights),
309  theHF1Grid(HF1Grid),
310  theHF1Weights(HF1Weights),
311  theHF2Grid(HF2Grid),
312  theHF2Weights(HF2Weights),
313  theEBweight(EBweight),
314  theEEweight(EEweight),
315  theHBweight(HBweight),
316  theHESweight(HESweight),
317  theHEDweight(HEDweight),
318  theHOweight(HOweight),
319  theHF1weight(HF1weight),
320  theHF2weight(HF2weight),
321  theEcutTower(EcutTower),
322  theEBSumThreshold(EBSumThreshold),
323  theEESumThreshold(EESumThreshold),
324  theEBEScale(50.),
325  theEEEScale(50.),
326  theHBEScale(50.),
327  theHESEScale(50.),
328  theHEDEScale(50.),
329  theHOEScale(50.),
330  theHF1EScale(50.),
331  theHF2EScale(50.),
332  theHcalTopology(nullptr),
333  theGeometry(nullptr),
334  theTowerConstituentsMap(nullptr),
335  theHcalAcceptSeverityLevel(0),
336  theRecoveredHcalHitsAreUsed(false),
337  theRecoveredEcalHitsAreUsed(false),
338  useRejectedHitsOnly(false),
339  theHcalAcceptSeverityLevelForRejectedHit(0),
340  useRejectedRecoveredHcalHits(0),
341  useRejectedRecoveredEcalHits(0),
343  theHOIsUsed(useHO),
344  // (momentum reconstruction algorithm)
345  theMomConstrMethod(momConstrMethod),
346  theMomHBDepth(momHBDepth),
347  theMomHEDepth(momHEDepth),
348  theMomEBDepth(momEBDepth),
349  theMomEEDepth(momEEDepth),
350  theHcalPhase(hcalPhase) {
351  // static int N = 0;
352  // std::cout << "VI Algo " << ++N << std::endl;
353  // nalgo=N;
354 }
355 
357 
359  const CaloTowerConstituentsMap* ctmap,
360  const HcalTopology* htopo,
361  const CaloGeometry* geo) {
362  theTowerTopology = cttopo;
363  theTowerConstituentsMap = ctmap;
364  theHcalTopology = htopo;
365  theGeometry = geo;
367 
368  //initialize ecal bad channel map
370 }
371 
373  theTowerMap.clear();
374  theTowerMapSize = 0;
375  //hcalDropChMap.clear();
376 }
377 
379  for (HBHERecHitCollection::const_iterator hbheItr = hbhe.begin(); hbheItr != hbhe.end(); ++hbheItr)
380  assignHitHcal(&(*hbheItr));
381 }
382 
384  for (HORecHitCollection::const_iterator hoItr = ho.begin(); hoItr != ho.end(); ++hoItr)
385  assignHitHcal(&(*hoItr));
386 }
387 
389  for (HFRecHitCollection::const_iterator hfItr = hf.begin(); hfItr != hf.end(); ++hfItr)
390  assignHitHcal(&(*hfItr));
391 }
392 
394  for (EcalRecHitCollection::const_iterator ecItr = ec.begin(); ecItr != ec.end(); ++ecItr)
395  assignHitEcal(&(*ecItr));
396 }
397 
398 // this method should not be used any more as the towers in the changed format
399 // can not be properly rescaled with the "rescale" method.
400 // "rescale was replaced by "rescaleTowers"
401 //
403  for (CaloTowerCollection::const_iterator ctcItr = ctc.begin(); ctcItr != ctc.end(); ++ctcItr) {
404  rescale(&(*ctcItr));
405  }
406 }
407 
409  // now copy this map into the final collection
410  result.reserve(theTowerMapSize);
411  // auto k=0U;
412  // if (!theEbHandle.isValid()) std::cout << "VI ebHandle not valid" << std::endl;
413  // if (!theEeHandle.isValid()) std::cout << "VI eeHandle not valid" << std::endl;
414 
415  for (auto const& mt : theTowerMap) {
416  // Convert only if there is at least one constituent in the metatower.
417  // The check of constituents size in the coverted tower is still needed!
418  if (!mt.empty()) {
419  convert(mt.id, mt, result);
420  } // ++k;}
421  }
422 
423  // assert(k==theTowerMapSize);
424  // std::cout << "VI TowerMap " << theTowerMapSize << " " << k << std::endl;
425 
426  theTowerMap.clear(); // save the memory
427  theTowerMapSize = 0;
428 }
429 
431  for (CaloTowerCollection::const_iterator ctcItr = ctc.begin(); ctcItr != ctc.end(); ++ctcItr) {
432  CaloTowerDetId twrId = ctcItr->id();
433  double newE_em = ctcItr->emEnergy();
434  double newE_had = ctcItr->hadEnergy();
435  double newE_outer = ctcItr->outerEnergy();
436 
437  double threshold = 0.0; // not used: we do not change thresholds
438  double weight = 1.0;
439 
440  // HF
441  if (ctcItr->ietaAbs() >= theTowerTopology->firstHFRing()) {
442  double E_short = 0.5 * newE_had; // from the definitions for HF
443  double E_long = newE_em + 0.5 * newE_had; //
444  // scale
445  E_long *= theHF1weight;
446  E_short *= theHF2weight;
447  // convert
448  newE_em = E_long - E_short;
449  newE_had = 2.0 * E_short;
450  }
451 
452  else { // barrel/endcap
453 
454  // find if its in EB, or EE; determine from first ecal constituent found
455  for (unsigned int iConst = 0; iConst < ctcItr->constituentsSize(); ++iConst) {
456  DetId constId = ctcItr->constituent(iConst);
457  if (constId.det() != DetId::Ecal)
458  continue;
460  newE_em *= weight;
461  break;
462  }
463  // HO
464  for (unsigned int iConst = 0; iConst < ctcItr->constituentsSize(); ++iConst) {
465  DetId constId = ctcItr->constituent(iConst);
466  if (constId.det() != DetId::Hcal)
467  continue;
468  if (HcalDetId(constId).subdet() != HcalOuter)
469  continue;
471  newE_outer *= weight;
472  break;
473  }
474  // HB/HE
475  for (unsigned int iConst = 0; iConst < ctcItr->constituentsSize(); ++iConst) {
476  DetId constId = ctcItr->constituent(iConst);
477  if (constId.det() != DetId::Hcal)
478  continue;
479  if (HcalDetId(constId).subdet() == HcalOuter)
480  continue;
482  newE_had *= weight;
483  if (ctcItr->ietaAbs() > theTowerTopology->firstHERing())
484  newE_outer *= weight;
485  break;
486  }
487 
488  } // barrel/endcap region
489 
490  // now make the new tower
491 
492  double newE_hadTot =
493  (theHOIsUsed && twrId.ietaAbs() <= theTowerTopology->lastHORing()) ? newE_had + newE_outer : newE_had;
494 
495  GlobalPoint emPoint = ctcItr->emPosition();
496  GlobalPoint hadPoint = ctcItr->emPosition();
497 
498  double f_em = 1.0 / cosh(emPoint.eta());
499  double f_had = 1.0 / cosh(hadPoint.eta());
500 
502 
503  if (ctcItr->ietaAbs() < theTowerTopology->firstHFRing()) {
504  if (newE_em > 0)
505  towerP4 += CaloTower::PolarLorentzVector(newE_em * f_em, emPoint.eta(), emPoint.phi(), 0);
506  if (newE_hadTot > 0)
507  towerP4 += CaloTower::PolarLorentzVector(newE_hadTot * f_had, hadPoint.eta(), hadPoint.phi(), 0);
508  } else {
509  double newE_tot = newE_em + newE_had;
510  // for HF we use common point for ecal, hcal shower positions regardless of the method
511  if (newE_tot > 0)
512  towerP4 += CaloTower::PolarLorentzVector(newE_tot * f_had, hadPoint.eta(), hadPoint.phi(), 0);
513  }
514 
515  CaloTower rescaledTower(twrId, newE_em, newE_had, newE_outer, -1, -1, towerP4, emPoint, hadPoint);
516  // copy the timings, have to convert back to int, 1 unit = 0.01 ns
517  rescaledTower.setEcalTime(int(ctcItr->ecalTime() * 100.0 + 0.5));
518  rescaledTower.setHcalTime(int(ctcItr->hcalTime() * 100.0 + 0.5));
519  //add topology info
520  rescaledTower.setHcalSubdet(theTowerTopology->lastHBRing(),
524 
525  std::vector<DetId> contains;
526  for (unsigned int iConst = 0; iConst < ctcItr->constituentsSize(); ++iConst) {
527  contains.push_back(ctcItr->constituent(iConst));
528  }
529  rescaledTower.addConstituents(contains);
530 
531  rescaledTower.setCaloTowerStatus(ctcItr->towerStatusWord());
532 
533  ctcResult.push_back(rescaledTower);
534 
535  } // end of loop over towers
536 }
537 
539  DetId detId = recHit->detid();
540  DetId detIdF(detId);
543 #ifdef EDM_ML_DEBUG
544  std::cout << "AssignHitHcal DetId " << HcalDetId(detId) << " Front " << HcalDetId(detIdF) << std::endl;
545 #endif
546  }
547 
548  unsigned int chStatusForCT = hcalChanStatusForCaloTower(recHit);
549 
550  // this is for skipping channls: mostly needed for the creation of
551  // bad towers from hits i the bad channel collections.
552  if (chStatusForCT == CaloTowersCreationAlgo::IgnoredChan)
553  return;
554 
555  double threshold, weight;
557 
558  double energy = recHit->energy(); // original RecHit energy is used to apply thresholds
559  double e = energy * weight; // energies scaled by user weight: used in energy assignments
560 
561  // SPECIAL handling of tower 28 merged depths --> half into tower 28 and half into tower 29
562  bool merge(false);
563  if (detIdF.det() == DetId::Hcal && HcalDetId(detIdF).subdet() == HcalEndcap &&
564  (theHcalPhase == 0 || theHcalPhase == 1) &&
565  //HcalDetId(detId).depth()==3 &&
566  HcalDetId(detIdF).ietaAbs() == theHcalTopology->lastHERing() - 1) {
568 #ifdef EDM_ML_DEBUG
569  std::cout << "Merge " << HcalDetId(detIdF) << ":" << merge << std::endl;
570 #endif
571  }
572  if (merge) {
574 
575  // bad channels are counted regardless of energy threshold
576 
577  if (chStatusForCT == CaloTowersCreationAlgo::BadChan) {
579  if (towerDetId.null())
580  return;
581  MetaTower& tower28 = find(towerDetId);
582  CaloTowerDetId towerDetId29(towerDetId.ieta() + towerDetId.zside(), towerDetId.iphi());
583  MetaTower& tower29 = find(towerDetId29);
584  tower28.numBadHcalCells += 1;
585  tower29.numBadHcalCells += 1;
586  }
587 
588  else if (0.5 * energy >= threshold) { // not bad channel: use energy if above threshold
589 
591  if (towerDetId.null())
592  return;
593  MetaTower& tower28 = find(towerDetId);
594  CaloTowerDetId towerDetId29(towerDetId.ieta() + towerDetId.zside(), towerDetId.iphi());
595  MetaTower& tower29 = find(towerDetId29);
596 
597  if (chStatusForCT == CaloTowersCreationAlgo::RecoveredChan) {
598  tower28.numRecHcalCells += 1;
599  tower29.numRecHcalCells += 1;
600  } else if (chStatusForCT == CaloTowersCreationAlgo::ProblematicChan) {
601  tower28.numProbHcalCells += 1;
602  tower29.numProbHcalCells += 1;
603  }
604 
605  // NOTE DIVIDE BY 2!!!
606  double e28 = 0.5 * e;
607  double e29 = 0.5 * e;
608 
609  tower28.E_had += e28;
610  tower28.E += e28;
611  std::pair<DetId, float> mc(detId, e28);
612  tower28.metaConstituents.push_back(mc);
613 
614  tower29.E_had += e29;
615  tower29.E += e29;
616  tower29.metaConstituents.push_back(mc);
617 
618  // time info: do not use in averaging if timing error is found: need
619  // full set of status info to implement: use only "good" channels for now
620 
621  if (chStatusForCT == CaloTowersCreationAlgo::GoodChan) {
622  tower28.hadSumTimeTimesE += (e28 * recHit->time());
623  tower28.hadSumEForTime += e28;
624  tower29.hadSumTimeTimesE += (e29 * recHit->time());
625  tower29.hadSumEForTime += e29;
626  }
627 
628  // store the energy in layer 3 also in E_outer
629  tower28.E_outer += e28;
630  tower29.E_outer += e29;
631  } // not a "bad" hit
632  } // end of special case
633 
634  else {
635  HcalDetId hcalDetId(detId);
636 
638 
639  if (hcalDetId.subdet() == HcalOuter) {
641  if (towerDetId.null())
642  return;
643  MetaTower& tower = find(towerDetId);
644 
645  if (chStatusForCT == CaloTowersCreationAlgo::BadChan) {
646  if (theHOIsUsed)
647  tower.numBadHcalCells += 1;
648  }
649 
650  else if (energy >= threshold) {
651  tower.E_outer += e; // store HO energy even if HO is not used
652  // add energy of the tower and/or flag if theHOIsUsed
653  if (theHOIsUsed) {
654  tower.E += e;
655 
656  if (chStatusForCT == CaloTowersCreationAlgo::RecoveredChan) {
657  tower.numRecHcalCells += 1;
658  } else if (chStatusForCT == CaloTowersCreationAlgo::ProblematicChan) {
659  tower.numProbHcalCells += 1;
660  }
661  } // HO is used
662 
663  // add HO to constituents even if it is not used: JetMET wants to keep these towers
664  std::pair<DetId, float> mc(detId, e);
665  tower.metaConstituents.push_back(mc);
666 
667  } // not a bad channel, energy above threshold
668 
669  } // HO hit
670 
671  // HF calculates EM fraction differently
672  else if (hcalDetId.subdet() == HcalForward) {
673  if (chStatusForCT == CaloTowersCreationAlgo::BadChan) {
675  if (towerDetId.null())
676  return;
677  MetaTower& tower = find(towerDetId);
678  tower.numBadHcalCells += 1;
679  }
680 
681  else if (energy >= threshold) {
683  if (towerDetId.null())
684  return;
685  MetaTower& tower = find(towerDetId);
686 
687  if (hcalDetId.depth() == 1) {
688  // long fiber, so E_EM = E(Long) - E(Short)
689  tower.E_em += e;
690  } else {
691  // short fiber, EHAD = 2 * E(Short)
692  tower.E_em -= e;
693  tower.E_had += 2. * e;
694  }
695  tower.E += e;
696  if (chStatusForCT == CaloTowersCreationAlgo::RecoveredChan) {
697  tower.numRecHcalCells += 1;
698  } else if (chStatusForCT == CaloTowersCreationAlgo::ProblematicChan) {
699  tower.numProbHcalCells += 1;
700  }
701 
702  // put the timing in HCAL -> have to check timing errors when available
703  // for now use only good channels
704  if (chStatusForCT == CaloTowersCreationAlgo::GoodChan) {
705  tower.hadSumTimeTimesE += (e * recHit->time());
706  tower.hadSumEForTime += e;
707  }
708 
709  std::pair<DetId, float> mc(detId, e);
710  tower.metaConstituents.push_back(mc);
711 
712  } // not a bad HF channel, energy above threshold
713 
714  } // HF hit
715 
716  else {
717  // HCAL situation normal in HB/HE
718  if (chStatusForCT == CaloTowersCreationAlgo::BadChan) {
720  if (towerDetId.null())
721  return;
722  MetaTower& tower = find(towerDetId);
723  tower.numBadHcalCells += 1;
724  } else if (energy >= threshold) {
726  if (towerDetId.null())
727  return;
728  MetaTower& tower = find(towerDetId);
729  tower.E_had += e;
730  tower.E += e;
731  if (chStatusForCT == CaloTowersCreationAlgo::RecoveredChan) {
732  tower.numRecHcalCells += 1;
733  } else if (chStatusForCT == CaloTowersCreationAlgo::ProblematicChan) {
734  tower.numProbHcalCells += 1;
735  }
736 
737  // Timing information: need specific accessors
738  // for now use only good channels
739  if (chStatusForCT == CaloTowersCreationAlgo::GoodChan) {
740  tower.hadSumTimeTimesE += (e * recHit->time());
741  tower.hadSumEForTime += e;
742  }
743  // store energy in highest depth for towers 18-27 (for electron,photon ID in endcap)
744  // also, store energy in HE part of tower 16 (for JetMET cleanup)
745  HcalDetId hcalDetId(detId);
746  if (hcalDetId.subdet() == HcalEndcap) {
747  if (theHcalPhase == 0) {
748  if ((hcalDetId.depth() == 2 && hcalDetId.ietaAbs() >= 18 && hcalDetId.ietaAbs() < 27) ||
749  (hcalDetId.depth() == 3 && hcalDetId.ietaAbs() == 27) ||
750  (hcalDetId.depth() == 3 && hcalDetId.ietaAbs() == 16)) {
751  tower.E_outer += e;
752  }
753  }
754  //combine depths in phase0-like way
755  else if (theHcalPhase == 1) {
756  if ((hcalDetId.depth() >= 3 && hcalDetId.ietaAbs() >= 18 && hcalDetId.ietaAbs() < 26) ||
757  (hcalDetId.depth() >= 4 && (hcalDetId.ietaAbs() == 26 || hcalDetId.ietaAbs() == 27)) ||
758  (hcalDetId.depth() == 3 && hcalDetId.ietaAbs() == 17) ||
759  (hcalDetId.depth() == 4 && hcalDetId.ietaAbs() == 16)) {
760  tower.E_outer += e;
761  }
762  }
763  }
764 
765  std::pair<DetId, float> mc(detId, e);
766  tower.metaConstituents.push_back(mc);
767 
768  } // not a "bad" channel, energy above threshold
769 
770  } // channel in HBHE (excluding twrs 28,29)
771 
772  } // recHit normal case (not in HE towers 28,29)
773 
774 } // end of assignHitHcal method
775 
777  DetId detId = recHit->detid();
778 
779  unsigned int chStatusForCT;
780  bool ecalIsBad = false;
781  std::tie(chStatusForCT, ecalIsBad) = ecalChanStatusForCaloTower(recHit);
782 
783  // this is for skipping channls: mostly needed for the creation of
784  // bad towers from hits i the bad channel collections.
785  if (chStatusForCT == CaloTowersCreationAlgo::IgnoredChan)
786  return;
787 
788  double threshold, weight;
790 
791  double energy = recHit->energy(); // original RecHit energy is used to apply thresholds
792  double e = energy * weight; // energies scaled by user weight: used in energy assignments
793 
795 
796  // For ECAL we count all bad channels after the metatower is complete
797 
798  // Include options for symmetric thresholds and cut on Et
799  // for ECAL RecHits
800 
801  bool passEmThreshold = false;
802 
803  if (detId.subdetId() == EcalBarrel) {
805  energy /= cosh((theGeometry->getGeometry(detId)->getPosition()).eta());
807  passEmThreshold = (fabs(energy) >= threshold);
808  else
809  passEmThreshold = (energy >= threshold);
810 
811  } else if (detId.subdetId() == EcalEndcap) {
813  energy /= cosh((theGeometry->getGeometry(detId)->getPosition()).eta());
815  passEmThreshold = (fabs(energy) >= threshold);
816  else
817  passEmThreshold = (energy >= threshold);
818  }
819 
821  if (towerDetId.null())
822  return;
823  MetaTower& tower = find(towerDetId);
824 
825  // count bad cells and avoid double counting with those from DB (Recovered are counted bad)
826 
827  // somehow misses some
828  // if ( (chStatusForCT == CaloTowersCreationAlgo::BadChan) & (!ecalIsBad) ) ++tower.numBadEcalCells;
829 
830  // a bit slower...
831  if (chStatusForCT == CaloTowersCreationAlgo::BadChan) {
832  auto thisEcalSevLvl = theEcalSevLvlAlgo->severityLevel(detId);
833  // check if the Ecal severity is ok to keep
834  auto sevit = std::find(theEcalSeveritiesToBeExcluded.begin(), theEcalSeveritiesToBeExcluded.end(), thisEcalSevLvl);
835  if (sevit == theEcalSeveritiesToBeExcluded.end())
836  ++tower.numBadEcalCells; // notinDB
837  }
838 
839  // if (chStatusForCT != CaloTowersCreationAlgo::BadChan && energy >= threshold) {
840  if (chStatusForCT != CaloTowersCreationAlgo::BadChan && passEmThreshold) {
841  tower.E_em += e;
842  tower.E += e;
843 
844  if (chStatusForCT == CaloTowersCreationAlgo::RecoveredChan) {
845  tower.numRecEcalCells += 1;
846  } else if (chStatusForCT == CaloTowersCreationAlgo::ProblematicChan) {
847  tower.numProbEcalCells += 1;
848  }
849 
850  // change when full status info is available
851  // for now use only good channels
852 
853  // add e>0 check (new options allow e<0)
854  if (chStatusForCT == CaloTowersCreationAlgo::GoodChan && e > 0) {
855  tower.emSumTimeTimesE += (e * recHit->time());
856  tower.emSumEForTime += e; // see above
857  }
858 
859  std::pair<DetId, float> mc(detId, e);
860  tower.metaConstituents.push_back(mc);
861  }
862 } // end of assignHitEcal method
863 
864 // This method is not flexible enough for the new CaloTower format.
865 // For now make a quick compatibility "fix" : WILL NOT WORK CORRECTLY with anything
866 // except the default simple p4 assignment!!!
867 // Must be rewritten for full functionality.
869  double threshold, weight;
870  CaloTowerDetId towerDetId = theTowerConstituentsMap->towerOf(ct->id());
871  if (towerDetId.null())
872  return;
873  MetaTower& tower = find(towerDetId);
874 
875  tower.E_em = 0.;
876  tower.E_had = 0.;
877  tower.E_outer = 0.;
878  for (unsigned int i = 0; i < ct->constituentsSize(); i++) {
879  DetId detId = ct->constituent(i);
881  DetId::Detector det = detId.det();
882  if (det == DetId::Ecal) {
883  tower.E_em = ct->emEnergy() * weight;
884  } else {
885  HcalDetId hcalDetId(detId);
886  if (hcalDetId.subdet() == HcalForward) {
887  if (hcalDetId.depth() == 1)
888  tower.E_em = ct->emEnergy() * weight;
889  if (hcalDetId.depth() == 2)
890  tower.E_had = ct->hadEnergy() * weight;
891  } else if (hcalDetId.subdet() == HcalOuter) {
892  tower.E_outer = ct->outerEnergy() * weight;
893  } else {
894  tower.E_had = ct->hadEnergy() * weight;
895  }
896  }
897  tower.E = tower.E_had + tower.E_em + tower.E_outer;
898 
899  // this is to be compliant with the new MetaTower setup
900  // used only for the default simple vector assignment
901  std::pair<DetId, float> mc(detId, 0);
902  tower.metaConstituents.push_back(mc);
903  }
904 
905  // preserve time inforamtion
906  tower.emSumTimeTimesE = ct->ecalTime();
907  tower.hadSumTimeTimesE = ct->hcalTime();
908  tower.emSumEForTime = 1.0;
909  tower.hadSumEForTime = 1.0;
910 }
911 
913  if (theTowerMap.empty()) {
915  }
916 
918 
919  if (mt.empty()) {
920  mt.id = detId;
921  mt.metaConstituents.reserve(detId.ietaAbs() < theTowerTopology->firstHFRing() ? 12 : 2);
922  ++theTowerMapSize;
923  }
924 
925  return mt;
926 }
927 
929  assert(id.rawId() != 0);
930 
931  double ecalThres = (id.ietaAbs() <= 17) ? (theEBSumThreshold) : (theEESumThreshold);
932  double E = mt.E;
933  double E_em = mt.E_em;
934  double E_had = mt.E_had;
935  double E_outer = mt.E_outer;
936 
937  // Note: E_outer is used to save HO energy OR energy in the outermost depths in endcap region
938  // In the methods with separate treatment of EM and HAD components:
939  // - HO is not used to determine direction, however HO energy is added to get "total had energy"
940  // => Check if the tower is within HO coverage before adding E_outer to the "total had" energy
941  // else the energy will be double counted
942  // When summing up the energy of the tower these checks are performed in the loops over RecHits
943 
944  std::vector<std::pair<DetId, float> > metaContains = mt.metaConstituents;
945  if (id.ietaAbs() < theTowerTopology->firstHFRing() && E_em < ecalThres) { // ignore EM threshold in HF
946  E -= E_em;
947  E_em = 0;
948  std::vector<std::pair<DetId, float> > metaContains_noecal;
949 
950  for (std::vector<std::pair<DetId, float> >::iterator i = metaContains.begin(); i != metaContains.end(); ++i)
951  if (i->first.det() != DetId::Ecal)
952  metaContains_noecal.push_back(*i);
953  metaContains.swap(metaContains_noecal);
954  }
955  if (id.ietaAbs() < theTowerTopology->firstHFRing() && E_had < theHcalThreshold) {
956  E -= E_had;
957 
959  E -= E_outer; // not subtracted before, think it should be done
960 
961  E_had = 0;
962  E_outer = 0;
963  std::vector<std::pair<DetId, float> > metaContains_nohcal;
964 
965  for (std::vector<std::pair<DetId, float> >::iterator i = metaContains.begin(); i != metaContains.end(); ++i)
966  if (i->first.det() != DetId::Hcal)
967  metaContains_nohcal.push_back(*i);
968  metaContains.swap(metaContains_nohcal);
969  }
970 
971  if (metaContains.empty())
972  return;
973 
974  if (missingHcalRescaleFactorForEcal > 0 && E_had == 0 && E_em > 0) {
975  auto match = hcalDropChMap.find(id);
976  if (match != hcalDropChMap.end() && match->second.second) {
977  E_had = missingHcalRescaleFactorForEcal * E_em;
978  E += E_had;
979  }
980  }
981 
982  double E_had_tot = (theHOIsUsed && id.ietaAbs() <= theTowerTopology->lastHORing()) ? E_had + E_outer : E_had;
983 
984  // create CaloTower using the selected algorithm
985 
986  GlobalPoint emPoint, hadPoint;
987 
988  // this is actually a 4D vector
989  Basic3DVectorF towerP4;
990  bool massless = true;
991  // float mass1=0;
992  float mass2 = 0;
993 
994  // conditional assignment of depths for barrel/endcap
995  // Some additional tuning may be required in the transitional region
996  // 14<|iEta|<19
997  double momEmDepth = 0.;
998  double momHadDepth = 0.;
999  if (id.ietaAbs() <= 17) {
1000  momHadDepth = theMomHBDepth;
1001  momEmDepth = theMomEBDepth;
1002  } else {
1003  momHadDepth = theMomHEDepth;
1004  momEmDepth = theMomEEDepth;
1005  }
1006 
1007  switch (theMomConstrMethod) {
1008  // FIXME : move to simple cartesian algebra
1009  case 0: { // Simple 4-momentum assignment
1010  GlobalPoint p = theTowerGeometry->getGeometry(id)->getPosition();
1011  towerP4 = p.basicVector().unit();
1012  towerP4[3] = 1.f; // energy
1013  towerP4 *= E;
1014 
1015  // double pf=1.0/cosh(p.eta());
1016  // if (E>0) towerP4 = CaloTower::PolarLorentzVector(E*pf, p.eta(), p.phi(), 0);
1017 
1018  emPoint = p;
1019  hadPoint = p;
1020  } // end case 0
1021  break;
1022 
1023  case 1: { // separate 4-vectors for ECAL, HCAL, add to get the 4-vector of the tower (=>tower has mass!)
1024  if (id.ietaAbs() < theTowerTopology->firstHFRing()) {
1025  Basic3DVectorF emP4;
1026  if (E_em > 0) {
1027  emPoint = emShwrPos(metaContains, momEmDepth, E_em);
1028  emP4 = emPoint.basicVector().unit();
1029  emP4[3] = 1.f; // energy
1030  towerP4 = emP4 * E_em;
1031 
1032  // double emPf = 1.0/cosh(emPoint.eta());
1033  // towerP4 += CaloTower::PolarLorentzVector(E_em*emPf, emPoint.eta(), emPoint.phi(), 0);
1034  }
1035  if ((E_had + E_outer) > 0) {
1036  massless = (E_em <= 0);
1037  hadPoint = hadShwrPos(id, momHadDepth);
1038  auto lP4 = hadPoint.basicVector().unit();
1039  lP4[3] = 1.f; // energy
1040  if (!massless) {
1041  auto diff = lP4 - emP4;
1042  mass2 = std::sqrt(E_em * E_had_tot * diff.mag2());
1043  }
1044  lP4 *= E_had_tot;
1045  towerP4 += lP4;
1046  /*
1047  if (!massless) {
1048  auto p = towerP4;
1049  double m2 = double(p[3]*p[3]) - double(p[0]*p[0])+double(p[1]*p[1])+double(p[2]*p[2]); mass1 = m2>0 ? std::sqrt(m2) : 0;
1050  }
1051  */
1052  // double hadPf = 1.0/cosh(hadPoint.eta());
1053  // if (E_had_tot>0) {
1054  // towerP4 += CaloTower::PolarLorentzVector(E_had_tot*hadPf, hadPoint.eta(), hadPoint.phi(), 0);
1055  // }
1056  }
1057  } else { // forward detector: use the CaloTower position
1058  GlobalPoint p = theTowerGeometry->getGeometry(id)->getPosition();
1059  towerP4 = p.basicVector().unit();
1060  towerP4[3] = 1.f; // energy
1061  towerP4 *= E;
1062  // double pf=1.0/cosh(p.eta());
1063  // if (E>0) towerP4 = CaloTower::PolarLorentzVector(E*pf, p.eta(), p.phi(), 0); // simple momentum assignment, same position
1064  emPoint = p;
1065  hadPoint = p;
1066  }
1067  } // end case 1
1068  break;
1069 
1070  case 2: { // use ECAL position for the tower (when E_cal>0), else default CaloTower position (massless tower)
1071  if (id.ietaAbs() < theTowerTopology->firstHFRing()) {
1072  if (E_em > 0)
1073  emPoint = emShwrLogWeightPos(metaContains, momEmDepth, E_em);
1074  else
1075  emPoint = theTowerGeometry->getGeometry(id)->getPosition();
1076  towerP4 = emPoint.basicVector().unit();
1077  towerP4[3] = 1.f; // energy
1078  towerP4 *= E;
1079 
1080  // double sumPf = 1.0/cosh(emPoint.eta());
1082 
1083  hadPoint = emPoint;
1084  } else { // forward detector: use the CaloTower position
1085  GlobalPoint p = theTowerGeometry->getGeometry(id)->getPosition();
1086  towerP4 = p.basicVector().unit();
1087  towerP4[3] = 1.f; // energy
1088  towerP4 *= E;
1089 
1090  // double pf=1.0/cosh(p.eta());
1091  // if (E>0) towerP4 = CaloTower::PolarLorentzVector(E*pf, p.eta(), p.phi(), 0); // simple momentum assignment, same position
1092  emPoint = p;
1093  hadPoint = p;
1094  }
1095  } // end case 2
1096  break;
1097 
1098  } // end of decision on p4 reconstruction method
1099 
1100  // insert in collection (remove and return if below threshold)
1101  if UNLIKELY ((towerP4[3] == 0) & (E_outer > 0)) {
1102  float val = theHOIsUsed ? 0 : 1E-9; // to keep backwards compatibility for theHOIsUsed == true
1103  collection.emplace_back(id,
1104  E_em,
1105  E_had,
1106  E_outer,
1107  -1,
1108  -1,
1109  CaloTower::PolarLorentzVector(val, hadPoint.eta(), hadPoint.phi(), 0),
1110  emPoint,
1111  hadPoint);
1112  } else {
1113  collection.emplace_back(
1114  id, E_em, E_had, E_outer, -1, -1, GlobalVector(towerP4), towerP4[3], mass2, emPoint, hadPoint);
1115  }
1116  auto& caloTower = collection.back();
1117 
1118  // if (!massless) std::cout << "massive " << id <<' ' << mass1 <<' ' << mass2 <<' ' << caloTower.mass() << std::endl;
1119  // std::cout << "CaloTowerVI " <<theMomConstrMethod <<' ' << id <<' '<< E_em <<' '<< E_had <<' '<< E_outer <<' '<< GlobalVector(towerP4) <<' '<< towerP4[3] <<' '<< emPoint <<' '<< hadPoint << std::endl;
1120  //if (towerP4[3]==0) std::cout << "CaloTowerVIzero " << theEcutTower << ' ' << collection.back().eta() <<' '<< collection.back().phi() << std::endl;
1121 
1122  if (caloTower.energy() < theEcutTower) {
1123  collection.pop_back();
1124  return;
1125  }
1126 
1127  // set the timings
1128  float ecalTime = (mt.emSumEForTime > 0) ? mt.emSumTimeTimesE / mt.emSumEForTime : -9999;
1129  float hcalTime = (mt.hadSumEForTime > 0) ? mt.hadSumTimeTimesE / mt.hadSumEForTime : -9999;
1130  caloTower.setEcalTime(compactTime(ecalTime));
1131  caloTower.setHcalTime(compactTime(hcalTime));
1132  //add topology info
1133  caloTower.setHcalSubdet(theTowerTopology->lastHBRing(),
1137 
1138  // set the CaloTower status word =====================================
1139  // Channels must be counter exclusively in the defined cathegories
1140  // "Bad" channels (not used in energy assignment) can be flagged during
1141  // CaloTower creation only if specified in the configuration file
1142 
1143  unsigned int numBadHcalChan = mt.numBadHcalCells;
1144  // unsigned int numBadEcalChan = mt.numBadEcalCells;
1145  unsigned int numBadEcalChan = 0; //
1146 
1147  unsigned int numRecHcalChan = mt.numRecHcalCells;
1148  unsigned int numRecEcalChan = mt.numRecEcalCells;
1149  unsigned int numProbHcalChan = mt.numProbHcalCells;
1150  unsigned int numProbEcalChan = mt.numProbEcalCells;
1151 
1152  // now add dead/off/... channels not used in RecHit reconstruction for HCAL
1153  HcalDropChMap::iterator dropChItr = hcalDropChMap.find(id);
1154  if (dropChItr != hcalDropChMap.end())
1155  numBadHcalChan += dropChItr->second.first;
1156 
1157  // for ECAL the number of all bad channels is obtained here -----------------------
1158 
1159  /*
1160  // old hyper slow algorithm
1161  // get all possible constituents of the tower
1162  std::vector<DetId> allConstituents = theTowerConstituentsMap->constituentsOf(id);
1163 
1164  for (std::vector<DetId>::iterator ac_it=allConstituents.begin();
1165  ac_it!=allConstituents.end(); ++ac_it) {
1166 
1167  if (ac_it->det()!=DetId::Ecal) continue;
1168 
1169  int thisEcalSevLvl = -999;
1170 
1171  if (ac_it->subdetId() == EcalBarrel && theEbHandle.isValid()) {
1172  thisEcalSevLvl = theEcalSevLvlAlgo->severityLevel( *ac_it, *theEbHandle);//, *theEcalChStatus);
1173  }
1174  else if (ac_it->subdetId() == EcalEndcap && theEeHandle.isValid()) {
1175  thisEcalSevLvl = theEcalSevLvlAlgo->severityLevel( *ac_it, *theEeHandle);//, *theEcalChStatus);
1176  }
1177 
1178  // check if the Ecal severity is ok to keep
1179  std::vector<int>::const_iterator sevit = std::find(theEcalSeveritiesToBeExcluded.begin(),
1180  theEcalSeveritiesToBeExcluded.end(),
1181  thisEcalSevLvl);
1182  if (sevit!=theEcalSeveritiesToBeExcluded.end()) {
1183  ++numBadEcalChan;
1184  }
1185 
1186  }
1187 
1188  // compare with fast version
1189 
1190  // hcal:
1191  int inEcals[2] = {0,0};
1192  for (std::vector<std::pair<DetId,float> >::iterator i=metaContains.begin(); i!=metaContains.end(); ++i) {
1193  DetId detId = i->first;
1194  if(detId.det() == DetId::Ecal){
1195  if( detId.subdetId()==EcalBarrel ) inEcals[0] =1;
1196  else if( detId.subdetId()==EcalEndcap ) inEcals[1] =1;
1197  }
1198  }
1199 
1200  auto numBadEcalChanNew = ecalBadChs[theTowerTopology->denseIndex(id)]+mt.numBadEcalCells; // - mt.numRecEcalCells
1201  if (int(numBadEcalChanNew)!=int(numBadEcalChan)) {
1202  std::cout << "VI wrong " << ((inEcals[1]==1) ? "EE" : "" ) << id << " " << numBadEcalChanNew << " " << numBadEcalChan
1203  << " " << mt.numBadEcalCells << " " << mt.numRecEcalCells << std::endl;
1204  }
1205  */
1206 
1207  numBadEcalChan = ecalBadChs[theTowerTopology->denseIndex(id)] + mt.numBadEcalCells; // - mt.numRecEcalCells
1208 
1209  //--------------------------------------------------------------------------------------
1210 
1211  caloTower.setCaloTowerStatus(
1212  numBadHcalChan, numBadEcalChan, numRecHcalChan, numRecEcalChan, numProbHcalChan, numProbEcalChan);
1213 
1214  double maxCellE = -999.0; // for storing the hottest cell E in the calotower
1215 
1216  std::vector<DetId> contains;
1217  contains.reserve(metaContains.size());
1218  for (std::vector<std::pair<DetId, float> >::iterator i = metaContains.begin(); i != metaContains.end(); ++i) {
1219  contains.push_back(i->first);
1220 
1221  if (maxCellE < i->second) {
1222  // need an extra check because of the funny towers that are empty except for the presence of an HO
1223  // hit in the constituents (JetMET wanted them saved)
1224  // This constituent is only used for storing the tower, but should not be concidered as a hot cell canditate for
1225  // configurations with useHO = false
1226 
1227  if (i->first.det() == DetId::Ecal) { // ECAL
1228  maxCellE = i->second;
1229  } else { // HCAL
1230  if (HcalDetId(i->first).subdet() != HcalOuter)
1231  maxCellE = i->second;
1232  else if (theHOIsUsed)
1233  maxCellE = i->second;
1234  }
1235 
1236  } // found higher E cell
1237 
1238  } // loop over matacontains
1239 
1240  caloTower.setConstituents(std::move(contains));
1241  caloTower.setHottestCellE(maxCellE);
1242 
1243  // std::cout << "CaloTowerVI " << nalgo << ' ' << caloTower.id() << ((inEcals[1]==1) ? "EE " : " " ) << caloTower.pt() << ' ' << caloTower.et() << ' ' << caloTower.mass() << ' '
1244  // << caloTower.constituentsSize() <<' '<< caloTower.towerStatusWord() << std::endl;
1245 }
1246 
1248  DetId::Detector det = detId.det();
1249  weight = 0; // in case the hit is not identified
1250 
1251  if (det == DetId::Ecal) {
1252  // may or may not be EB. We'll find out.
1253 
1254  EcalSubdetector subdet = (EcalSubdetector)(detId.subdetId());
1255  if (subdet == EcalBarrel) {
1257  weight = theEBweight;
1258  if (weight <= 0.) {
1259  ROOT::Math::Interpolator my(theEBGrid, theEBWeights, ROOT::Math::Interpolation::kAKIMA);
1260  weight = my.Eval(theEBEScale);
1261  }
1262  } else if (subdet == EcalEndcap) {
1264  weight = theEEweight;
1265  if (weight <= 0.) {
1266  ROOT::Math::Interpolator my(theEEGrid, theEEWeights, ROOT::Math::Interpolation::kAKIMA);
1267  weight = my.Eval(theEEEScale);
1268  }
1269  }
1270  } else if (det == DetId::Hcal) {
1271  HcalDetId hcalDetId(detId);
1272  HcalSubdetector subdet = hcalDetId.subdet();
1273  int depth = hcalDetId.depth();
1274 
1275  if (subdet == HcalBarrel) {
1276  if (hcalCuts == nullptr) { // this means cutsFromDB is false
1278  } else { // hcalCuts is not nullptr, i.e. cutsFromDB is true
1279  const HcalPFCut* item = hcalCuts->getValues(hcalDetId.rawId());
1280  threshold = item->noiseThreshold();
1281  }
1282  weight = theHBweight;
1283  if (weight <= 0.) {
1284  ROOT::Math::Interpolator my(theHBGrid, theHBWeights, ROOT::Math::Interpolation::kAKIMA);
1285  weight = my.Eval(theHBEScale);
1286  }
1287  }
1288 
1289  else if (subdet == HcalEndcap) {
1290  // check if it's single or double tower
1291  if (hcalDetId.ietaAbs() < theHcalTopology->firstHEDoublePhiRing()) {
1292  if (hcalCuts == nullptr) { // this means cutsFromDB is false
1294  } else { // hcalCuts is not nullptr, i.e. cutsFromDB is true
1295  const HcalPFCut* item = hcalCuts->getValues(hcalDetId.rawId());
1296  threshold = item->noiseThreshold();
1297  }
1298  weight = theHESweight;
1299  if (weight <= 0.) {
1300  ROOT::Math::Interpolator my(theHESGrid, theHESWeights, ROOT::Math::Interpolation::kAKIMA);
1301  weight = my.Eval(theHESEScale);
1302  }
1303  } else {
1304  if (hcalCuts == nullptr) { // this means cutsFromDB is false
1306  } else { // hcalCuts is not nullptr, i.e. cutsFromDB is true
1307  const HcalPFCut* item = hcalCuts->getValues(hcalDetId.rawId());
1308  threshold = item->noiseThreshold();
1309  }
1310  weight = theHEDweight;
1311  if (weight <= 0.) {
1312  ROOT::Math::Interpolator my(theHEDGrid, theHEDWeights, ROOT::Math::Interpolation::kAKIMA);
1313  weight = my.Eval(theHEDEScale);
1314  }
1315  }
1316  }
1317 
1318  else if (subdet == HcalOuter) {
1319  //check if it's ring 0 or +1 or +2 or -1 or -2
1320  if (hcalDetId.ietaAbs() <= 4)
1322  else if (hcalDetId.ieta() < 0) {
1323  // set threshold for ring -1 or -2
1324  threshold = (hcalDetId.ietaAbs() <= 10) ? theHOthresholdMinus1 : theHOthresholdMinus2;
1325  } else {
1326  // set threshold for ring +1 or +2
1327  threshold = (hcalDetId.ietaAbs() <= 10) ? theHOthresholdPlus1 : theHOthresholdPlus2;
1328  }
1329  weight = theHOweight;
1330  if (weight <= 0.) {
1331  ROOT::Math::Interpolator my(theHOGrid, theHOWeights, ROOT::Math::Interpolation::kAKIMA);
1332  weight = my.Eval(theHOEScale);
1333  }
1334  }
1335 
1336  else if (subdet == HcalForward) {
1337  if (hcalDetId.depth() == 1) {
1339  weight = theHF1weight;
1340  if (weight <= 0.) {
1341  ROOT::Math::Interpolator my(theHF1Grid, theHF1Weights, ROOT::Math::Interpolation::kAKIMA);
1342  weight = my.Eval(theHF1EScale);
1343  }
1344  } else {
1346  weight = theHF2weight;
1347  if (weight <= 0.) {
1348  ROOT::Math::Interpolator my(theHF2Grid, theHF2Weights, ROOT::Math::Interpolation::kAKIMA);
1349  weight = my.Eval(theHF2EScale);
1350  }
1351  }
1352  }
1353  } else {
1354  edm::LogError("CaloTowersCreationAlgo") << "Bad cell: " << det << std::endl;
1355  }
1356 }
1357 
1359  if (scale > 0.00001)
1360  *&theEBEScale = scale;
1361  else
1362  *&theEBEScale = 50.;
1363 }
1364 
1366  if (scale > 0.00001)
1367  *&theEEEScale = scale;
1368  else
1369  *&theEEEScale = 50.;
1370 }
1371 
1373  if (scale > 0.00001)
1374  *&theHBEScale = scale;
1375  else
1376  *&theHBEScale = 50.;
1377 }
1378 
1380  if (scale > 0.00001)
1381  *&theHESEScale = scale;
1382  else
1383  *&theHESEScale = 50.;
1384 }
1385 
1387  if (scale > 0.00001)
1388  *&theHEDEScale = scale;
1389  else
1390  *&theHEDEScale = 50.;
1391 }
1392 
1394  if (scale > 0.00001)
1395  *&theHOEScale = scale;
1396  else
1397  *&theHOEScale = 50.;
1398 }
1399 
1401  if (scale > 0.00001)
1402  *&theHF1EScale = scale;
1403  else
1404  *&theHF1EScale = 50.;
1405 }
1406 
1408  if (scale > 0.00001)
1409  *&theHF2EScale = scale;
1410  else
1411  *&theHF2EScale = 50.;
1412 }
1413 
1415  auto cellGeometry = theGeometry->getGeometry(detId);
1416  GlobalPoint point = cellGeometry->getPosition(); // face of the cell
1417 
1418  if (fracDepth <= 0)
1419  return point;
1420  if (fracDepth > 1)
1421  fracDepth = 1;
1422 
1423  const GlobalPoint& backPoint = cellGeometry->getBackPoint();
1424  point += fracDepth * (backPoint - point);
1425 
1426  return point;
1427 }
1428 
1430  // same code as above
1431  return emCrystalShwrPos(detId, fracDepth);
1432 }
1433 
1434 GlobalPoint CaloTowersCreationAlgo::hadShwrPos(const std::vector<std::pair<DetId, float> >& metaContains,
1435  float fracDepth,
1436  double hadE) {
1437  // this is based on available RecHits, can lead to different actual depths if
1438  // hits in multi-depth towers are not all there
1439 #ifdef EDM_ML_DEBUG
1440  std::cout << "hadShwrPos called with " << metaContains.size() << " elements and energy " << hadE << ":" << fracDepth
1441  << std::endl;
1442 #endif
1443  if (hadE <= 0)
1444  return GlobalPoint(0, 0, 0);
1445 
1446  double hadX = 0.0;
1447  double hadY = 0.0;
1448  double hadZ = 0.0;
1449 
1450  int nConst = 0;
1451 
1452  std::vector<std::pair<DetId, float> >::const_iterator mc_it = metaContains.begin();
1453  for (; mc_it != metaContains.end(); ++mc_it) {
1454  if (mc_it->first.det() != DetId::Hcal)
1455  continue;
1456  // do not use HO for deirection calculations for now
1457  if (HcalDetId(mc_it->first).subdet() == HcalOuter)
1458  continue;
1459  ++nConst;
1460 
1461  GlobalPoint p = hadSegmentShwrPos(mc_it->first, fracDepth);
1462 
1463  // longitudinal segmentation: do not weight by energy,
1464  // get the geometrical position
1465  hadX += p.x();
1466  hadY += p.y();
1467  hadZ += p.z();
1468  }
1469 
1470  return GlobalPoint(hadX / nConst, hadY / nConst, hadZ / nConst);
1471 }
1472 
1474  // set depth using geometry of cells that are associated with the
1475  // tower (regardless if they have non-zero energies)
1476 
1477  // if (hadE <= 0) return GlobalPoint(0, 0, 0);
1478 #ifdef EDM_ML_DEBUG
1479  std::cout << "hadShwrPos " << towerId << " frac " << fracDepth << std::endl;
1480 #endif
1481  if (fracDepth < 0)
1482  fracDepth = 0;
1483  else if (fracDepth > 1)
1484  fracDepth = 1;
1485 
1486  GlobalPoint point(0, 0, 0);
1487 
1488  int iEta = towerId.ieta();
1489  int iPhi = towerId.iphi();
1490 
1491  HcalDetId frontCellId, backCellId;
1492 
1493  if (towerId.ietaAbs() >= theTowerTopology->firstHFRing()) {
1494  // forward, take the geometry for long fibers
1495  frontCellId = HcalDetId(HcalForward, towerId.zside() * theTowerTopology->convertCTtoHcal(abs(iEta)), iPhi, 1);
1496  backCellId = HcalDetId(HcalForward, towerId.zside() * theTowerTopology->convertCTtoHcal(abs(iEta)), iPhi, 1);
1497  } else {
1498  //use constituents map
1499  std::vector<DetId> items = theTowerConstituentsMap->constituentsOf(towerId);
1500  int frontDepth = 1000;
1501  int backDepth = -1000;
1502  for (unsigned i = 0; i < items.size(); i++) {
1503  if (items[i].det() != DetId::Hcal)
1504  continue;
1505  HcalDetId hid(items[i]);
1506  if (hid.subdet() == HcalOuter)
1507  continue;
1508  if (!theHcalTopology->validHcal(hid, 2))
1509  continue;
1510 
1511  if (theHcalTopology->idFront(hid).depth() < frontDepth) {
1512  frontCellId = hid;
1513  frontDepth = theHcalTopology->idFront(hid).depth();
1514  }
1515  if (theHcalTopology->idBack(hid).depth() > backDepth) {
1516  backCellId = hid;
1517  backDepth = theHcalTopology->idBack(hid).depth();
1518  }
1519  }
1520 #ifdef EDM_ML_DEBUG
1521  std::cout << "Front " << frontCellId << " Back " << backCellId << " Depths " << frontDepth << ":" << backDepth
1522  << std::endl;
1523 #endif
1524  //fix for tower 28/29 - no tower 29 at highest depths
1525  if (towerId.ietaAbs() == theTowerTopology->lastHERing() && (theHcalPhase == 0 || theHcalPhase == 1)) {
1526  CaloTowerDetId towerId28(towerId.ieta() - towerId.zside(), towerId.iphi());
1527  std::vector<DetId> items28 = theTowerConstituentsMap->constituentsOf(towerId28);
1528 #ifdef EDM_ML_DEBUG
1529  std::cout << towerId28 << " with " << items28.size() << " constituents:";
1530  for (unsigned k = 0; k < items28.size(); ++k)
1531  if (items28[k].det() == DetId::Hcal)
1532  std::cout << " " << HcalDetId(items28[k]);
1533  std::cout << std::endl;
1534 #endif
1535  for (unsigned i = 0; i < items28.size(); i++) {
1536  if (items28[i].det() != DetId::Hcal)
1537  continue;
1538  HcalDetId hid(items28[i]);
1539  if (hid.subdet() == HcalOuter)
1540  continue;
1541 
1542  if (theHcalTopology->idBack(hid).depth() > backDepth) {
1543  backCellId = hid;
1544  backDepth = theHcalTopology->idBack(hid).depth();
1545  }
1546  }
1547  }
1548 #ifdef EDM_ML_DEBUG
1549  std::cout << "Back " << backDepth << " ID " << backCellId << std::endl;
1550 #endif
1551  }
1552  point = hadShwPosFromCells(DetId(frontCellId), DetId(backCellId), fracDepth);
1553 
1554  return point;
1555 }
1556 
1557 GlobalPoint CaloTowersCreationAlgo::hadShwPosFromCells(DetId frontCellId, DetId backCellId, float fracDepth) {
1558  // uses the "front" and "back" cells
1559  // to determine the axis. point set by the predefined depth.
1560 
1561  HcalDetId hid1(frontCellId), hid2(backCellId);
1563  hid1 = theHcalTopology->idFront(frontCellId);
1564 #ifdef EDM_ML_DEBUG
1565  std::cout << "Front " << HcalDetId(frontCellId) << " " << hid1 << "\n";
1566 #endif
1567  hid2 = theHcalTopology->idBack(backCellId);
1568 #ifdef EDM_ML_DEBUG
1569  std::cout << "Back " << HcalDetId(backCellId) << " " << hid2 << "\n";
1570 #endif
1571  }
1572 
1573  auto frontCellGeometry = theGeometry->getGeometry(DetId(hid1));
1574  auto backCellGeometry = theGeometry->getGeometry(DetId(hid2));
1575 
1576  GlobalPoint point = frontCellGeometry->getPosition();
1577  const GlobalPoint& backPoint = backCellGeometry->getBackPoint();
1578 
1579  point += fracDepth * (backPoint - point);
1580 
1581  return point;
1582 }
1583 
1584 GlobalPoint CaloTowersCreationAlgo::emShwrPos(const std::vector<std::pair<DetId, float> >& metaContains,
1585  float fracDepth,
1586  double emE) {
1587  if (emE <= 0)
1588  return GlobalPoint(0, 0, 0);
1589 
1590  double emX = 0.0;
1591  double emY = 0.0;
1592  double emZ = 0.0;
1593 
1594  double eSum = 0;
1595 
1596  std::vector<std::pair<DetId, float> >::const_iterator mc_it = metaContains.begin();
1597  for (; mc_it != metaContains.end(); ++mc_it) {
1598  if (mc_it->first.det() != DetId::Ecal)
1599  continue;
1600  GlobalPoint p = emCrystalShwrPos(mc_it->first, fracDepth);
1601  double e = mc_it->second;
1602 
1603  if (e > 0) {
1604  emX += p.x() * e;
1605  emY += p.y() * e;
1606  emZ += p.z() * e;
1607  eSum += e;
1608  }
1609  }
1610 
1611  return GlobalPoint(emX / eSum, emY / eSum, emZ / eSum);
1612 }
1613 
1614 GlobalPoint CaloTowersCreationAlgo::emShwrLogWeightPos(const std::vector<std::pair<DetId, float> >& metaContains,
1615  float fracDepth,
1616  double emE) {
1617  double emX = 0.0;
1618  double emY = 0.0;
1619  double emZ = 0.0;
1620 
1621  double weight = 0;
1622  double sumWeights = 0;
1623  double sumEmE = 0; // add crystals with E/E_EM > 1.5%
1624  double crystalThresh = 0.015 * emE;
1625 
1626  std::vector<std::pair<DetId, float> >::const_iterator mc_it = metaContains.begin();
1627  for (; mc_it != metaContains.end(); ++mc_it) {
1628  if (mc_it->second < 0)
1629  continue;
1630  if (mc_it->first.det() == DetId::Ecal && mc_it->second > crystalThresh)
1631  sumEmE += mc_it->second;
1632  }
1633 
1634  for (mc_it = metaContains.begin(); mc_it != metaContains.end(); ++mc_it) {
1635  if (mc_it->first.det() != DetId::Ecal || mc_it->second < crystalThresh)
1636  continue;
1637 
1638  GlobalPoint p = emCrystalShwrPos(mc_it->first, fracDepth);
1639 
1640  weight = 4.2 + log(mc_it->second / sumEmE);
1641  sumWeights += weight;
1642 
1643  emX += p.x() * weight;
1644  emY += p.y() * weight;
1645  emZ += p.z() * weight;
1646  }
1647 
1648  return GlobalPoint(emX / sumWeights, emY / sumWeights, emZ / sumWeights);
1649 }
1650 
1652  const float timeUnit = 0.01; // discretization (ns)
1653 
1654  if (time > 300.0)
1655  return 30000;
1656  if (time < -300.0)
1657  return -30000;
1658 
1659  return int(time / timeUnit + 0.5);
1660 }
1661 
1662 //========================================================
1663 //
1664 // Bad/anomolous cell handling
1665 
1667  // This method fills the map of number of dead channels for the calotower,
1668  // The key of the map is CaloTowerDetId.
1669  // By definition these channels are not going to be in the RecHit collections.
1670  hcalDropChMap.clear();
1671  std::vector<DetId> allChanInStatusCont = theHcalChStatus->getAllChannels();
1672 
1673 #ifdef EDM_ML_DEBUG
1674  std::cout << "DropChMap with " << allChanInStatusCont.size() << " channels" << std::endl;
1675 #endif
1676  for (std::vector<DetId>::iterator it = allChanInStatusCont.begin(); it != allChanInStatusCont.end(); ++it) {
1677  const uint32_t dbStatusFlag = theHcalChStatus->getValues(*it)->getValue();
1678  if (theHcalSevLvlComputer->dropChannel(dbStatusFlag)) {
1680 
1682 
1683  hcalDropChMap[twrId].first += 1;
1684 
1685  HcalDetId hid(*it);
1686 
1687  // special case for tower 29: if HCAL hit is in depth 3 add to twr 29 as well
1688  if (hid.subdet() == HcalEndcap && (theHcalPhase == 0 || theHcalPhase == 1) &&
1689  hid.ietaAbs() == theHcalTopology->lastHERing() - 1) {
1690  bool merge = theHcalTopology->mergedDepth29(hid);
1691  if (merge) {
1692  CaloTowerDetId twrId29(twrId.ieta() + twrId.zside(), twrId.iphi());
1693  hcalDropChMap[twrId29].first += 1;
1694  }
1695  }
1696  }
1697  }
1698  // now I know how many bad channels, but I also need to know if there's any good ones
1700  for (auto& pair : hcalDropChMap) {
1701  if (pair.second.first == 0)
1702  continue; // unexpected, but just in case
1703  int ngood = 0, nbad = 0;
1704  for (DetId id : theTowerConstituentsMap->constituentsOf(pair.first)) {
1705  if (id.det() != DetId::Hcal)
1706  continue;
1707  HcalDetId hid(id);
1708  if (hid.subdet() != HcalBarrel && hid.subdet() != HcalEndcap)
1709  continue;
1710  const uint32_t dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
1711  if (dbStatusFlag == 0 || !theHcalSevLvlComputer->dropChannel(dbStatusFlag)) {
1712  ngood += 1;
1713  } else {
1714  nbad += 1; // recount, since pair.second.first may include HO
1715  }
1716  }
1717  if (nbad > 0 && nbad >= ngood) {
1718  //uncomment for debug (may be useful to tune the criteria above)
1719  //CaloTowerDetId id(pair.first);
1720  //std::cout << "CaloTower at ieta = " << id.ieta() << ", iphi " << id.iphi() << ": set Hcal as not efficient (ngood =" << ngood << ", nbad = " << nbad << ")" << std::endl;
1721  pair.second.second = true;
1722  }
1723  }
1724  }
1725 }
1726 
1728  // std::cout << "VI making EcalBadChs ";
1729 
1730  // for ECAL the number of all bad channels is obtained here -----------------------
1731 
1732  for (auto ind = 0U; ind < theTowerTopology->sizeForDenseIndexing(); ++ind) {
1733  auto& numBadEcalChan = ecalBadChs[ind];
1734  numBadEcalChan = 0;
1735  auto id = theTowerTopology->detIdFromDenseIndex(ind);
1736 
1737  // this is utterly slow... (can be optmized if really needed)
1738 
1739  // get all possible constituents of the tower
1740  std::vector<DetId> allConstituents = theTowerConstituentsMap->constituentsOf(id);
1741 
1742  for (std::vector<DetId>::iterator ac_it = allConstituents.begin(); ac_it != allConstituents.end(); ++ac_it) {
1743  if (ac_it->det() != DetId::Ecal)
1744  continue;
1745 
1746  auto thisEcalSevLvl = theEcalSevLvlAlgo->severityLevel(*ac_it);
1747 
1748  // check if the Ecal severity is ok to keep
1749  std::vector<int>::const_iterator sevit =
1751  if (sevit != theEcalSeveritiesToBeExcluded.end()) {
1752  ++numBadEcalChan;
1753  }
1754  }
1755 
1756  // if (0!=numBadEcalChan) std::cout << id << ":" << numBadEcalChan << ", ";
1757  }
1758 
1759  /*
1760  int tot=0;
1761  for (auto ind=0U; ind<theTowerTopology->sizeForDenseIndexing(); ++ind) {
1762  if (ecalBadChs[ind]!=0) ++tot;
1763  }
1764  std::cout << " | " << tot << std::endl;
1765  */
1766 }
1767 
1769 
1771  HcalDetId hid(hit->detid());
1772  DetId id = theHcalTopology->idFront(hid);
1773 #ifdef EDM_ML_DEBUG
1774  std::cout << "ChanStatusForCaloTower for " << hid << " to " << HcalDetId(id) << std::endl;
1775 #endif
1776  const uint32_t recHitFlag = hit->flags();
1777  const uint32_t dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
1778 
1779  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
1780  bool isRecovered = theHcalSevLvlComputer->recoveredRecHit(id, recHitFlag);
1781 
1782  // For use with hits rejected in the default reconstruction
1783  if (useRejectedHitsOnly) {
1784  if (!isRecovered) {
1788  // this hit was either already accepted or is worse than
1789  } else {
1791  // skip recovered hits either because they were already used or because there was an explicit instruction
1793  } else if (useRejectedRecoveredHcalHits) {
1795  }
1796 
1797  } // recovered channels
1798 
1799  // clasify channels as problematic: no good hits are supposed to be present in the
1800  // extra rechit collections
1802 
1803  } // treatment of rejected hits
1804 
1805  // this is for the regular reconstruction sequence
1806 
1807  if (severityLevel == 0)
1809 
1810  if (isRecovered) {
1812  } else {
1815  } else {
1817  }
1818  }
1819 }
1820 
1822  // const DetId id = hit->detid();
1823 
1824  // uint16_t dbStatus = theEcalChStatus->find(id)->getStatusCode();
1825  // uint32_t rhFlags = hit->flags();
1826  // int severityLevel = theEcalSevLvlAlgo->severityLevel(rhFlags, dbStatus);
1827  // The methods above will become private and cannot be usef for flagging ecal spikes.
1828  // Use the recommended interface - we leave the parameters for spilke removal to be specified by ECAL.
1829 
1830  // int severityLevel = 999;
1831 
1832  EcalRecHit const& rh = *reinterpret_cast<EcalRecHit const*>(hit);
1834 
1835  // if (id.subdetId() == EcalBarrel) severityLevel = theEcalSevLvlAlgo->severityLevel( id, *theEbHandle);//, *theEcalChStatus);
1836  // else if (id.subdetId() == EcalEndcap) severityLevel = theEcalSevLvlAlgo->severityLevel( id, *theEeHandle);//, *theEcalChStatus);
1837 
1838  // there should be no other ECAL types used in this reconstruction
1839 
1840  // The definition of ECAL severity levels uses categories that
1841  // are similar to the defined for CaloTower. (However, the categorization
1842  // for CaloTowers depends on the specified maximum acceptabel severity and therefore cannnot
1843  // be exact correspondence between the two. ECAL has additional categories describing modes of failure.)
1844  // This approach is different from the initial idea and from
1845  // the implementation for HCAL. Still make the logic similar to HCAL so that one has the ability to
1846  // exclude problematic channels as defined by ECAL.
1847  // For definitions of ECAL severity levels see RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h
1848 
1849  bool isBad = (severityLevel == EcalSeverityLevel::kBad);
1850 
1851  bool isRecovered = (severityLevel == EcalSeverityLevel::kRecovered);
1852 
1853  // check if the severity is compatible with our configuration
1854  // This applies to the "default" tower cleaning
1855  std::vector<int>::const_iterator sevit =
1857  bool accepted = (sevit == theEcalSeveritiesToBeExcluded.end());
1858 
1859  // For use with hits that were rejected in the regular reconstruction:
1860  // This is for creating calotowers with lower level of cleaning by merging
1861  // the information from the default towers and a collection of towers created from
1862  // bad rechits
1863 
1864  if (useRejectedHitsOnly) {
1865  if (!isRecovered) {
1866  if (accepted || std::find(theEcalSeveritiesToBeUsedInBadTowers.begin(),
1869  return std::make_tuple(CaloTowersCreationAlgo::IgnoredChan, isBad);
1870  // this hit was either already accepted, or is not eligible for inclusion
1871  } else {
1873  // skip recovered hits either because they were already used or because there was an explicit instruction
1874  return std::make_tuple(CaloTowersCreationAlgo::IgnoredChan, isBad);
1875  ;
1876  } else if (useRejectedRecoveredEcalHits) {
1877  return std::make_tuple(CaloTowersCreationAlgo::RecoveredChan, isBad);
1878  }
1879 
1880  } // recovered channels
1881 
1882  // clasify channels as problematic
1883  return std::make_tuple(CaloTowersCreationAlgo::ProblematicChan, isBad);
1884 
1885  } // treatment of rejected hits
1886 
1887  // for normal reconstruction
1889  return std::make_tuple(CaloTowersCreationAlgo::GoodChan, false);
1890 
1891  if (isRecovered) {
1892  return std::make_tuple(
1894  } else {
1895  return std::make_tuple(accepted ? CaloTowersCreationAlgo::ProblematicChan : CaloTowersCreationAlgo::BadChan, isBad);
1896  }
1897 }
std::vector< double > theHBGrid
float hcalTime() const
Definition: CaloTower.h:197
bool mergedDepth29(HcalDetId id) const
Definition: HcalTopology.h:111
Definition: merge.py:1
bool getMergePositionFlag() const
Definition: HcalTopology.h:167
const HcalChannelQuality * theHcalChStatus
uint32_t sizeForDenseIndexing() const
std::vector< std::pair< DetId, float > > metaConstituents
bool contains(EventRange const &lh, EventID const &rh)
Definition: EventRange.cc:37
double outerEnergy() const
Definition: CaloTower.h:132
int ietaAbs() const
get the absolute value of the tower ieta
int merge(int argc, char *argv[])
Definition: DiMuonVmerge.cc:27
void getThresholdAndWeight(const DetId &detId, double &threshold, double &weight) const
helper method to look up the appropriate threshold & weight
constexpr int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
void setCaloTowerStatus(unsigned int numBadHcalChan, unsigned int numBadEcalChan, unsigned int numRecHcalChan, unsigned int numRecEcalChan, unsigned int numProbHcalChan, unsigned int numProbEcalChan)
Definition: CaloTower.cc:221
bool recoveredRecHit(const DetId &myid, const uint32_t &myflag) const
CaloTowerDetId towerOf(const DetId &id) const
Get the tower id for this det id (or null if not known)
std::vector< double > theHESGrid
std::vector< double > theHEDGrid
GlobalPoint emCrystalShwrPos(DetId detId, float fracDepth)
void assignHitHcal(const CaloRecHit *recHit)
std::vector< int > theEcalSeveritiesToBeUsedInBadTowers
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: LeafCandidate.h:25
T eta() const
Definition: PV3DBase.h:73
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< T >::const_iterator const_iterator
void push_back(T const &t)
Definition: weight.py:1
void finish(CaloTowerCollection &destCollection)
std::vector< double > theHOWeights
MetaTower & find(const CaloTowerDetId &id)
looks for a given tower in the internal cache. If it can&#39;t find it, it makes it.
void setThresFromDB(const HcalPFCuts *cuts)
Log< level::Error, false > LogError
std::vector< double > theEEGrid
void rescale(const CaloTower *ct)
int lastHFRing() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
void setHcalSubdet(int lastHB, int lastHE, int lastHF, int lastHO)
Definition: CaloTower.h:85
const Item * getValues(DetId fId, bool throwOnFail=true) const
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
std::vector< double > theEEWeights
int lastHERing() const
int lastHORing() const
std::vector< double > theHESWeights
U second(std::pair< T, U > const &p)
void setGeometry(const CaloTowerTopology *cttopo, const CaloTowerConstituentsMap *ctmap, const HcalTopology *htopo, const CaloGeometry *geo)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ietaAbs(uint32_t id)
constexpr bool null() const
is this a null id ?
Definition: DetId.h:59
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
const CaloSubdetectorGeometry * theTowerGeometry
void addConstituents(const std::vector< DetId > &ids)
Definition: CaloTower.cc:199
unsigned int hcalChanStatusForCaloTower(const CaloRecHit *hit)
CaloTowerDetId detIdFromDenseIndex(uint32_t din) const
int firstHERing() const
std::vector< double > theHF2Grid
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:170
std::vector< double > theHEDWeights
unsigned towerId(DetId const &, EcalElectronicsMapping const *)
void rescaleTowers(const CaloTowerCollection &ctInput, CaloTowerCollection &ctResult)
static const int SubdetId
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
T sqrt(T t)
Definition: SSEVec.h:19
const CaloGeometry * theGeometry
int convertCTtoHcal(int ct_ieta) const
HcalDetId idBack(const HcalDetId &id) const
Definition: HcalTopology.h:171
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:60
std::vector< double > theHOGrid
void assignHitEcal(const EcalRecHit *recHit)
adds a single hit to the tower
std::vector< DetId > getAllChannels() const
const CaloTowerConstituentsMap * theTowerConstituentsMap
void convert(const CaloTowerDetId &id, const MetaTower &mt, CaloTowerCollection &collection)
const CaloTowerTopology * theTowerTopology
std::vector< double > theHF1Weights
const_iterator begin() const
GlobalPoint hadSegmentShwrPos(DetId detId, float fracDepth)
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
uint32_t getValue() const
bool theHOIsUsed
only affects energy and ET calculation. HO is still recorded in the tower
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
DetId constituent(size_t i) const
Definition: CaloTower.h:126
void process(const HBHERecHitCollection &hbhe)
const_iterator end() const
std::vector< double > theHBWeights
const HcalSeverityLevelComputer * theHcalSevLvlComputer
Definition: DetId.h:17
GlobalPoint hadShwrPos(const std::vector< std::pair< DetId, float >> &metaContains, float fracDepth, double hadE)
int iphi() const
get the tower iphi
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
double emEnergy() const
Definition: CaloTower.h:130
Detector
Definition: DetId.h:24
int firstHEDoublePhiRing() const
Definition: HcalTopology.h:101
GlobalPoint hadShwPosFromCells(DetId frontCell, DetId backCell, float fracDepth)
size_t constituentsSize() const
Definition: CaloTower.h:125
double hadEnergy() const
Definition: CaloTower.h:131
int ieta() const
get the tower ieta
const HcalTopology * theHcalTopology
bool validHcal(const HcalDetId &id) const
std::vector< int > theEcalSeveritiesToBeExcluded
std::vector< unsigned short > ecalBadChs
float ecalTime() const
Definition: CaloTower.h:196
void setEcalTime(int t)
Definition: CaloTower.h:83
HcalDetId mergedDepthDetId(const HcalDetId &id) const
Definition: HcalTopology.h:166
GlobalPoint emShwrPos(const std::vector< std::pair< DetId, float >> &metaContains, float fracDepth, double totEmE)
const EcalSeverityLevelAlgo * theEcalSevLvlAlgo
int zside() const
get the z-side of the tower (1/-1)
int lastHERing() const
Definition: HcalTopology.h:94
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
int firstHFRing() const
std::vector< double > theHF1Grid
uint32_t denseIndex(const DetId &id) const
std::vector< double > theEBGrid
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
Basic3DVector unit() const
std::vector< double > theHF2Weights
#define UNLIKELY(x)
Definition: Likely.h:21
unsigned int theHcalAcceptSeverityLevelForRejectedHit
std::tuple< unsigned int, bool > ecalChanStatusForCaloTower(const EcalRecHit *hit)
void setHcalTime(int t)
Definition: CaloTower.h:84
EcalSubdetector
std::vector< double > theEBWeights
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
CaloTowerDetId id() const
Definition: CaloTower.h:123
std::vector< DetId > constituentsOf(const CaloTowerDetId &id) const
Get the constituent detids for this tower id ( not yet implemented )
def move(src, dest)
Definition: eostools.py:511
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalPoint emShwrLogWeightPos(const std::vector< std::pair< DetId, float >> &metaContains, float fracDepth, double totEmE)
int lastHBRing() const
bool dropChannel(const uint32_t &mystatus) const
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164