CMS 3D CMS Logo

HcalDigitizer.cc
Go to the documentation of this file.
30 #include <boost/foreach.hpp>
31 #include <memory>
32 
33 //#define EDM_ML_DEBUG
34 
36  : conditionsToken_(iC.esConsumes()),
37  topoToken_(iC.esConsumes()),
38  topoZToken_(iC.esConsumes()),
39  hcalTimeSlew_delay_token_(iC.esConsumes(edm::ESInputTag("", "HBHE"))),
40  theGeometryToken(iC.esConsumes()),
41  theRecNumberToken(iC.esConsumes()),
42  qieTypesToken_(iC.esConsumes()),
43  theGeometry(nullptr),
44  theRecNumber(nullptr),
45  theParameterMap(ps),
46  theShapes(),
47  theHBHEResponse(std::make_unique<CaloHitResponse>(&theParameterMap, &theShapes)),
48  theHBHESiPMResponse(std::make_unique<HcalSiPMHitResponse>(
49  &theParameterMap, &theShapes, ps.getParameter<bool>("HcalPreMixStage1"), true)),
50  theHOResponse(std::make_unique<CaloHitResponse>(&theParameterMap, &theShapes)),
51  theHOSiPMResponse(std::make_unique<HcalSiPMHitResponse>(
52  &theParameterMap, &theShapes, ps.getParameter<bool>("HcalPreMixStage1"), false)),
53  theHFResponse(std::make_unique<CaloHitResponse>(&theParameterMap, &theShapes)),
54  theHFQIE10Response(std::make_unique<CaloHitResponse>(
55  &theParameterMap, &theShapes, ps.getParameter<bool>("HcalPreMixStage1"), true)),
56  theZDCResponse(std::make_unique<CaloHitResponse>(
57  &theParameterMap, &theShapes, ps.getParameter<bool>("HcalPreMixStage1"), false)),
58  theHBHEAmplifier(nullptr),
59  theHFAmplifier(nullptr),
60  theHOAmplifier(nullptr),
61  theZDCAmplifier(nullptr),
62  theHFQIE10Amplifier(nullptr),
63  theHBHEQIE11Amplifier(nullptr),
64  theIonFeedback(nullptr),
65  theCoderFactory(nullptr),
66  theHBHEElectronicsSim(nullptr),
67  theHFElectronicsSim(nullptr),
68  theHOElectronicsSim(nullptr),
69  theZDCElectronicsSim(nullptr),
70  theHFQIE10ElectronicsSim(nullptr),
71  theHBHEQIE11ElectronicsSim(nullptr),
72  theHBHEHitFilter(),
73  theHBHEQIE11HitFilter(),
74  theHFHitFilter(),
75  theHFQIE10HitFilter(),
76  theHOHitFilter(),
77  theHOSiPMHitFilter(),
78  theZDCHitFilter(),
79  theHBHEDigitizer(nullptr),
80  theHODigitizer(nullptr),
81  theHOSiPMDigitizer(nullptr),
82  theHFDigitizer(nullptr),
83  theZDCDigitizer(nullptr),
84  theHFQIE10Digitizer(nullptr),
85  theHBHEQIE11Digitizer(nullptr),
86  theRelabeller(nullptr),
87  isZDC(true),
88  isHCAL(true),
89  zdcgeo(true),
90  hbhegeo(true),
91  hogeo(true),
92  hfgeo(true),
93  doHFWindow_(ps.getParameter<bool>("doHFWindow")),
94  killHE_(ps.getParameter<bool>("killHE")),
95  debugCS_(ps.getParameter<bool>("debugCaloSamples")),
96  ignoreTime_(ps.getParameter<bool>("ignoreGeantTime")),
97  injectTestHits_(ps.getParameter<bool>("injectTestHits")),
98  hitsProducer_(ps.getParameter<std::string>("hitsProducer")),
99  theHOSiPMCode(ps.getParameter<edm::ParameterSet>("ho").getParameter<int>("siPMCode")),
100  deliveredLumi(0.),
101  agingFlagHB(ps.getParameter<bool>("HBDarkening")),
102  agingFlagHE(ps.getParameter<bool>("HEDarkening")),
103  zdcToken_(iC.consumes(edm::InputTag(hitsProducer_, "ZDCHITS"))),
104  hcalToken_(iC.consumes(edm::InputTag(hitsProducer_, "HcalHits"))),
105  m_HBDarkening(nullptr),
106  m_HEDarkening(nullptr),
107  m_HFRecalibration(nullptr),
108  injectedHitsEnergy_(ps.getParameter<std::vector<double>>("injectTestHitsEnergy")),
109  injectedHitsTime_(ps.getParameter<std::vector<double>>("injectTestHitsTime")),
110  injectedHitsCells_(ps.getParameter<std::vector<int>>("injectTestHitsCells")) {
111  if (agingFlagHB) {
113  }
114  if (agingFlagHE) {
116  }
117  if (theHOSiPMCode == 2) {
118  mcParamsToken_ = iC.esConsumes();
119  }
120 
121  bool doNoise = ps.getParameter<bool>("doNoise");
122 
123  bool PreMix1 = ps.getParameter<bool>("HcalPreMixStage1"); // special threshold/pedestal treatment
124  bool PreMix2 = ps.getParameter<bool>("HcalPreMixStage2"); // special threshold/pedestal treatment
125  bool doEmpty = ps.getParameter<bool>("doEmpty");
126  deliveredLumi = ps.getParameter<double>("DelivLuminosity");
127  bool agingFlagHF = ps.getParameter<bool>("HFDarkening");
128  double minFCToDelay = ps.getParameter<double>("minFCToDelay");
129 
130  if (PreMix1 && PreMix2) {
131  throw cms::Exception("Configuration") << "HcalDigitizer cannot operate in PreMixing digitization and "
132  "PreMixing\n"
133  "digi combination modes at the same time. Please set one mode to "
134  "False\n"
135  "in the configuration file.";
136  }
137 
138  // need to make copies, because they might get different noise generators
139  theHBHEAmplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
140  theHFAmplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
141  theHOAmplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
142  theZDCAmplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
143  theHFQIE10Amplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
144  theHBHEQIE11Amplifier = std::make_unique<HcalAmplifier>(&theParameterMap, doNoise, PreMix1, PreMix2);
145 
146  theCoderFactory = std::make_unique<HcalCoderFactory>(HcalCoderFactory::DB);
147 
149  std::make_unique<HcalElectronicsSim>(&theParameterMap, theHBHEAmplifier.get(), theCoderFactory.get(), PreMix1);
151  std::make_unique<HcalElectronicsSim>(&theParameterMap, theHFAmplifier.get(), theCoderFactory.get(), PreMix1);
153  std::make_unique<HcalElectronicsSim>(&theParameterMap, theHOAmplifier.get(), theCoderFactory.get(), PreMix1);
155  std::make_unique<HcalElectronicsSim>(&theParameterMap, theZDCAmplifier.get(), theCoderFactory.get(), PreMix1);
157  std::make_unique<HcalElectronicsSim>(&theParameterMap,
158  theHFQIE10Amplifier.get(),
159  theCoderFactory.get(),
160  PreMix1); // should this use a different coder factory?
162  std::make_unique<HcalElectronicsSim>(&theParameterMap,
163  theHBHEQIE11Amplifier.get(),
164  theCoderFactory.get(),
165  PreMix1); // should this use a different coder factory?
166 
167  bool doHOHPD = (theHOSiPMCode != 1);
168  bool doHOSiPM = (theHOSiPMCode != 0);
169  if (doHOHPD) {
170  theHOResponse = std::make_unique<CaloHitResponse>(&theParameterMap, &theShapes);
171  theHOResponse->setHitFilter(&theHOHitFilter);
172  theHODigitizer = std::make_unique<HODigitizer>(theHOResponse.get(), theHOElectronicsSim.get(), doEmpty);
173  }
174  if (doHOSiPM) {
175  theHOSiPMResponse->setHitFilter(&theHOSiPMHitFilter);
176  theHOSiPMDigitizer = std::make_unique<HODigitizer>(theHOSiPMResponse.get(), theHOElectronicsSim.get(), doEmpty);
177  }
178 
179  theHBHEResponse->setHitFilter(&theHBHEHitFilter);
181 
182  // QIE8 and QIE11 can coexist in HBHE
184  std::make_unique<QIE11Digitizer>(theHBHESiPMResponse.get(), theHBHEQIE11ElectronicsSim.get(), doEmpty);
185  theHBHEDigitizer = std::make_unique<HBHEDigitizer>(theHBHEResponse.get(), theHBHEElectronicsSim.get(), doEmpty);
186 
187  bool doTimeSlew = ps.getParameter<bool>("doTimeSlew");
188  // initialize: they won't be called later if flag is set
189  hcalTimeSlew_delay_ = nullptr;
190  theTimeSlewSim.reset(nullptr);
191  if (doTimeSlew) {
192  // no time slewing for HF
193  theTimeSlewSim = std::make_unique<HcalTimeSlewSim>(&theParameterMap, minFCToDelay);
194  theHBHEAmplifier->setTimeSlewSim(theTimeSlewSim.get());
195  theHBHEQIE11Amplifier->setTimeSlewSim(theTimeSlewSim.get());
196  theHOAmplifier->setTimeSlewSim(theTimeSlewSim.get());
197  theZDCAmplifier->setTimeSlewSim(theTimeSlewSim.get());
198  }
199 
200  theHFResponse->setHitFilter(&theHFHitFilter);
201  theHFQIE10Response->setHitFilter(&theHFQIE10HitFilter);
202  theZDCResponse->setHitFilter(&theZDCHitFilter);
203 
204  // QIE8 and QIE10 can coexist in HF
206  std::make_unique<QIE10Digitizer>(theHFQIE10Response.get(), theHFQIE10ElectronicsSim.get(), doEmpty);
207  theHFDigitizer = std::make_unique<HFDigitizer>(theHFResponse.get(), theHFElectronicsSim.get(), doEmpty);
208 
209  // temporary move until Run3 ZDC will be added
210  theZDCDigitizer = std::make_unique<ZDCDigitizer>(theZDCResponse.get(), theZDCElectronicsSim.get(), doEmpty);
211 
212  testNumbering_ = ps.getParameter<bool>("TestNumbering");
213  // edm::LogVerbatim("HcalSim") << "Flag to see if Hit Relabeller to be initiated " << testNumbering_;
214  if (testNumbering_)
215  theRelabeller = std::make_unique<HcalHitRelabeller>(ps.getParameter<bool>("doNeutralDensityFilter"));
216 
217  if (ps.getParameter<bool>("doIonFeedback") && theHBHEResponse) {
218  theIonFeedback = std::make_unique<HPDIonFeedbackSim>(ps, &theShapes);
219  theHBHEResponse->setPECorrection(theIonFeedback.get());
220  if (ps.getParameter<bool>("doThermalNoise")) {
221  theHBHEAmplifier->setIonFeedbackSim(theIonFeedback.get());
222  }
223  }
224 
225  // option to save CaloSamples as event product for debugging
226  if (debugCS_) {
227  if (theHBHEDigitizer)
228  theHBHEDigitizer->setDebugCaloSamples(true);
230  theHBHEQIE11Digitizer->setDebugCaloSamples(true);
231  if (theHODigitizer)
232  theHODigitizer->setDebugCaloSamples(true);
233  if (theHOSiPMDigitizer)
234  theHOSiPMDigitizer->setDebugCaloSamples(true);
235  if (theHFDigitizer)
236  theHFDigitizer->setDebugCaloSamples(true);
238  theHFQIE10Digitizer->setDebugCaloSamples(true);
239  if (theZDCDigitizer)
240  theZDCDigitizer->setDebugCaloSamples(true);
241  }
242 
243  // option to ignore Geant time distribution in SimHits, for debugging
244  if (ignoreTime_) {
245  theHBHEResponse->setIgnoreGeantTime(ignoreTime_);
246  theHBHESiPMResponse->setIgnoreGeantTime(ignoreTime_);
247  theHOResponse->setIgnoreGeantTime(ignoreTime_);
248  theHOSiPMResponse->setIgnoreGeantTime(ignoreTime_);
249  theHFResponse->setIgnoreGeantTime(ignoreTime_);
250  theHFQIE10Response->setIgnoreGeantTime(ignoreTime_);
251  theZDCResponse->setIgnoreGeantTime(ignoreTime_);
252  }
253 
254  if (agingFlagHF)
255  m_HFRecalibration = std::make_unique<HFRecalibration>(ps.getParameter<edm::ParameterSet>("HFRecalParameterBlock"));
256 }
257 
259 
261  noiseGenerator->setParameterMap(&theParameterMap);
262  noiseGenerator->setElectronicsSim(theHBHEElectronicsSim.get());
263  if (theHBHEDigitizer)
264  theHBHEDigitizer->setNoiseSignalGenerator(noiseGenerator);
265  theHBHEAmplifier->setNoiseSignalGenerator(noiseGenerator);
266 }
267 
269  noiseGenerator->setParameterMap(&theParameterMap);
270  noiseGenerator->setElectronicsSim(theHBHEQIE11ElectronicsSim.get());
272  theHBHEQIE11Digitizer->setNoiseSignalGenerator(noiseGenerator);
273  theHBHEQIE11Amplifier->setNoiseSignalGenerator(noiseGenerator);
274 }
275 
277  noiseGenerator->setParameterMap(&theParameterMap);
278  noiseGenerator->setElectronicsSim(theHFElectronicsSim.get());
279  if (theHFDigitizer)
280  theHFDigitizer->setNoiseSignalGenerator(noiseGenerator);
281  theHFAmplifier->setNoiseSignalGenerator(noiseGenerator);
282 }
283 
285  noiseGenerator->setParameterMap(&theParameterMap);
286  noiseGenerator->setElectronicsSim(theHFQIE10ElectronicsSim.get());
288  theHFQIE10Digitizer->setNoiseSignalGenerator(noiseGenerator);
289  theHFQIE10Amplifier->setNoiseSignalGenerator(noiseGenerator);
290 }
291 
293  noiseGenerator->setParameterMap(&theParameterMap);
294  noiseGenerator->setElectronicsSim(theHOElectronicsSim.get());
295  if (theHODigitizer)
296  theHODigitizer->setNoiseSignalGenerator(noiseGenerator);
297  if (theHOSiPMDigitizer)
298  theHOSiPMDigitizer->setNoiseSignalGenerator(noiseGenerator);
299  theHOAmplifier->setNoiseSignalGenerator(noiseGenerator);
300 }
301 
303  noiseGenerator->setParameterMap(&theParameterMap);
304  noiseGenerator->setElectronicsSim(theZDCElectronicsSim.get());
305  if (theZDCDigitizer) {
306  theZDCDigitizer->setNoiseSignalGenerator(noiseGenerator);
307  theZDCAmplifier->setNoiseSignalGenerator(noiseGenerator);
308  }
309 }
310 
312  setup(eventSetup);
313 
314  // get the appropriate gains, noises, & widths for this event
316 
318 
319  theHBHEAmplifier->setDbService(conditions);
320  theHFAmplifier->setDbService(conditions);
321  theHOAmplifier->setDbService(conditions);
322  theZDCAmplifier->setDbService(conditions);
323  theHFQIE10Amplifier->setDbService(conditions);
324  theHBHEQIE11Amplifier->setDbService(conditions);
325 
326  theHFQIE10ElectronicsSim->setDbService(conditions);
327  theHBHEQIE11ElectronicsSim->setDbService(conditions);
328 
329  theCoderFactory->setDbService(conditions);
331 
332  // initialize hits
333  if (theHBHEDigitizer)
334  theHBHEDigitizer->initializeHits();
336  theHBHEQIE11Digitizer->initializeHits();
337  if (theHODigitizer)
338  theHODigitizer->initializeHits();
339  if (theHOSiPMDigitizer)
340  theHOSiPMDigitizer->initializeHits();
342  theHFQIE10Digitizer->initializeHits();
343  if (theHFDigitizer)
344  theHFDigitizer->initializeHits();
345  if (theZDCDigitizer)
346  theZDCDigitizer->initializeHits();
347 }
348 
349 void HcalDigitizer::accumulateCaloHits(edm::Handle<std::vector<PCaloHit>> const &hcalHandle,
350  edm::Handle<std::vector<PCaloHit>> const &zdcHandle,
351  int bunchCrossing,
352  CLHEP::HepRandomEngine *engine,
353  const HcalTopology *htopoP,
354  const ZdcTopology *ztopoP) {
355  // Step A: pass in inputs, and accumulate digis
356  if (isHCAL) {
357  std::vector<PCaloHit> hcalHitsOrig = *hcalHandle.product();
358  if (injectTestHits_)
359  hcalHitsOrig = injectedHits_;
360  std::vector<PCaloHit> hcalHits;
361  hcalHits.reserve(hcalHitsOrig.size());
362 
363  // evaluate darkening before relabeling
364  if (testNumbering_) {
366  darkening(hcalHitsOrig);
367  }
368  // Relabel PCaloHits if necessary
369  edm::LogInfo("HcalDigitizer") << "Calling Relabeller";
370  theRelabeller->process(hcalHitsOrig);
371  }
372 
373  // eliminate bad hits
374  for (unsigned int i = 0; i < hcalHitsOrig.size(); i++) {
375  DetId id(hcalHitsOrig[i].id());
376  HcalDetId hid(id);
377  if (!htopoP->validHcal(hid)) {
378  edm::LogError("HcalDigitizer") << "bad hcal id found in digitizer. Skipping " << id.rawId() << " " << hid;
379  continue;
380  } else if (hid.subdet() == HcalForward && !doHFWindow_ && hcalHitsOrig[i].depth() != 0) {
381  // skip HF window hits unless desired
382  continue;
383  } else if (killHE_ && hid.subdet() == HcalEndcap) {
384  // remove HE hits if asked for (phase 2)
385  continue;
386  } else {
387 #ifdef EDM_ML_DEBUG
388  edm::LogVerbatim("HcalSim") << "HcalDigitizer format " << hid.oldFormat() << " for " << hid;
389 #endif
390  DetId newid = DetId(hid.newForm());
391 #ifdef EDM_ML_DEBUG
392  edm::LogVerbatim("HcalSim") << "Hit " << i << " out of " << hcalHits.size() << " " << std::hex << id.rawId()
393  << " --> " << newid.rawId() << std::dec << " " << HcalDetId(newid.rawId()) << '\n';
394 #endif
395  hcalHitsOrig[i].setID(newid.rawId());
396  hcalHits.push_back(hcalHitsOrig[i]);
397  }
398  }
399 
400  if (hbhegeo) {
401  if (theHBHEDigitizer)
402  theHBHEDigitizer->add(hcalHits, bunchCrossing, engine);
404  theHBHEQIE11Digitizer->add(hcalHits, bunchCrossing, engine);
405  }
406 
407  if (hogeo) {
408  if (theHODigitizer)
409  theHODigitizer->add(hcalHits, bunchCrossing, engine);
410  if (theHOSiPMDigitizer)
411  theHOSiPMDigitizer->add(hcalHits, bunchCrossing, engine);
412  }
413 
414  if (hfgeo) {
415  if (theHFDigitizer)
416  theHFDigitizer->add(hcalHits, bunchCrossing, engine);
418  theHFQIE10Digitizer->add(hcalHits, bunchCrossing, engine);
419  }
420  } else {
421  edm::LogInfo("HcalDigitizer") << "We don't have HCAL hit collection available ";
422  }
423 
424  if (isZDC) {
425  if (zdcgeo && theZDCDigitizer) {
426  std::vector<PCaloHit> zdcHitsOrig = *zdcHandle.product();
427  std::vector<PCaloHit> zdcHits;
428  zdcHits.reserve(zdcHitsOrig.size());
429  // eliminate bad hits
430  for (unsigned int i = 0; i < zdcHitsOrig.size(); i++) {
431  DetId id(zdcHitsOrig[i].id());
432  HcalZDCDetId hid(id);
433  if (!ztopoP->valid(hid)) {
434  edm::LogError("HcalDigitizer") << "bad zdc id found in digitizer. Skipping " << std::hex << id.rawId()
435  << std::dec << " " << hid;
436  continue;
437  }
438  zdcHits.push_back(zdcHitsOrig[i]);
439 #ifdef EDM_ML_DEBUG
440  edm::LogVerbatim("HcalSim") << "Hit " << i << " out of " << zdcHitsOrig.size() << " " << std::hex << id.rawId()
441  << " " << hid;
442 #endif
443  }
444  theZDCDigitizer->add(zdcHits, bunchCrossing, engine);
445  }
446  } else {
447  edm::LogInfo("HcalDigitizer") << "We don't have ZDC hit collection available ";
448  }
449 }
450 
451 void HcalDigitizer::accumulate(edm::Event const &e, edm::EventSetup const &eventSetup, CLHEP::HepRandomEngine *engine) {
452  // Step A: Get Inputs
453  const edm::Handle<std::vector<PCaloHit>> &zdcHandle = e.getHandle(zdcToken_);
454  isZDC = zdcHandle.isValid();
455 
456  const edm::Handle<std::vector<PCaloHit>> &hcalHandle = e.getHandle(hcalToken_);
457  isHCAL = hcalHandle.isValid() or injectTestHits_;
458 
459  const HcalTopology *htopoP = &eventSetup.getData(topoToken_);
460  const ZdcTopology *ztopoP = &eventSetup.getData(topoZToken_);
461 
462  accumulateCaloHits(hcalHandle, zdcHandle, 0, engine, htopoP, ztopoP);
463 }
464 
467  CLHEP::HepRandomEngine *engine) {
468  // Step A: Get Inputs
469  edm::InputTag zdcTag(hitsProducer_, "ZDCHITS");
471  e.getByLabel(zdcTag, zdcHandle);
472  isZDC = zdcHandle.isValid();
473 
474  edm::InputTag hcalTag(hitsProducer_, "HcalHits");
476  e.getByLabel(hcalTag, hcalHandle);
477  isHCAL = hcalHandle.isValid();
478 
479  const HcalTopology *htopoP = &eventSetup.getData(topoToken_);
480  const ZdcTopology *ztopoP = &eventSetup.getData(topoZToken_);
481 
482  accumulateCaloHits(hcalHandle, zdcHandle, e.bunchCrossing(), engine, htopoP, ztopoP);
483 }
484 
485 void HcalDigitizer::finalizeEvent(edm::Event &e, const edm::EventSetup &eventSetup, CLHEP::HepRandomEngine *engine) {
486  // Step B: Create empty output
487  std::unique_ptr<HBHEDigiCollection> hbheResult(new HBHEDigiCollection());
488  std::unique_ptr<HODigiCollection> hoResult(new HODigiCollection());
489  std::unique_ptr<HFDigiCollection> hfResult(new HFDigiCollection());
490  std::unique_ptr<ZDCDigiCollection> zdcResult(new ZDCDigiCollection());
491  std::unique_ptr<QIE10DigiCollection> hfQIE10Result(new QIE10DigiCollection(
492  !theHFQIE10DetIds.empty() ? theHFQIE10Response.get()->getReadoutFrameSize(theHFQIE10DetIds[0])
494  std::unique_ptr<QIE11DigiCollection> hbheQIE11Result(new QIE11DigiCollection(
495  !theHBHEQIE11DetIds.empty() ? theHBHESiPMResponse.get()->getReadoutFrameSize(theHBHEQIE11DetIds[0]) :
496  // theParameterMap->simParameters(theHBHEQIE11DetIds[0]).readoutFrameSize()
497  // :
499 
500  // Step C: Invoke the algorithm, getting back outputs.
501  if (isHCAL && hbhegeo) {
502  if (theHBHEDigitizer)
503  theHBHEDigitizer->run(*hbheResult, engine);
505  theHBHEQIE11Digitizer->run(*hbheQIE11Result, engine);
506  }
507  if (isHCAL && hogeo) {
508  if (theHODigitizer)
509  theHODigitizer->run(*hoResult, engine);
510  if (theHOSiPMDigitizer)
511  theHOSiPMDigitizer->run(*hoResult, engine);
512  }
513  if (isHCAL && hfgeo) {
514  if (theHFDigitizer)
515  theHFDigitizer->run(*hfResult, engine);
517  theHFQIE10Digitizer->run(*hfQIE10Result, engine);
518  }
519  if (isZDC && zdcgeo && theZDCDigitizer) {
520  theZDCDigitizer->run(*zdcResult, engine);
521  }
522 
523  edm::LogInfo("HcalDigitizer") << "HCAL HBHE digis : " << hbheResult->size();
524  edm::LogInfo("HcalDigitizer") << "HCAL HO digis : " << hoResult->size();
525  edm::LogInfo("HcalDigitizer") << "HCAL HF digis : " << hfResult->size();
526  edm::LogInfo("HcalDigitizer") << "HCAL ZDC digis : " << zdcResult->size();
527  edm::LogInfo("HcalDigitizer") << "HCAL HF QIE10 digis : " << hfQIE10Result->size();
528  edm::LogInfo("HcalDigitizer") << "HCAL HBHE QIE11 digis : " << hbheQIE11Result->size();
529 
530 #ifdef EDM_ML_DEBUG
531  edm::LogVerbatim("HcalSim") << "\nHCAL HBHE digis : " << hbheResult->size();
532  edm::LogVerbatim("HcalSim") << "HCAL HO digis : " << hoResult->size();
533  edm::LogVerbatim("HcalSim") << "HCAL HF digis : " << hfResult->size();
534  edm::LogVerbatim("HcalSim") << "HCAL ZDC digis : " << zdcResult->size();
535  edm::LogVerbatim("HcalSim") << "HCAL HF QIE10 digis : " << hfQIE10Result->size();
536  edm::LogVerbatim("HcalSim") << "HCAL HBHE QIE11 digis : " << hbheQIE11Result->size();
537 #endif
538 
539  // Step D: Put outputs into event
540  e.put(std::move(hbheResult));
541  e.put(std::move(hoResult));
542  e.put(std::move(hfResult));
543  e.put(std::move(zdcResult));
544  e.put(std::move(hfQIE10Result), "HFQIE10DigiCollection");
545  e.put(std::move(hbheQIE11Result), "HBHEQIE11DigiCollection");
546 
547  if (debugCS_) {
548  std::unique_ptr<CaloSamplesCollection> csResult(new CaloSamplesCollection());
549  // smush together all the results
550  if (theHBHEDigitizer)
551  csResult->insert(
552  csResult->end(), theHBHEDigitizer->getCaloSamples().begin(), theHBHEDigitizer->getCaloSamples().end());
554  csResult->insert(csResult->end(),
555  theHBHEQIE11Digitizer->getCaloSamples().begin(),
556  theHBHEQIE11Digitizer->getCaloSamples().end());
557  if (theHODigitizer)
558  csResult->insert(
559  csResult->end(), theHODigitizer->getCaloSamples().begin(), theHODigitizer->getCaloSamples().end());
560  if (theHOSiPMDigitizer)
561  csResult->insert(
562  csResult->end(), theHOSiPMDigitizer->getCaloSamples().begin(), theHOSiPMDigitizer->getCaloSamples().end());
563  if (theHFDigitizer)
564  csResult->insert(
565  csResult->end(), theHFDigitizer->getCaloSamples().begin(), theHFDigitizer->getCaloSamples().end());
567  csResult->insert(
568  csResult->end(), theHFQIE10Digitizer->getCaloSamples().begin(), theHFQIE10Digitizer->getCaloSamples().end());
569  if (theZDCDigitizer)
570  csResult->insert(
571  csResult->end(), theZDCDigitizer->getCaloSamples().begin(), theZDCDigitizer->getCaloSamples().end());
572  e.put(std::move(csResult), "HcalSamples");
573  }
574 
575  if (injectTestHits_) {
576  std::unique_ptr<edm::PCaloHitContainer> pcResult(new edm::PCaloHitContainer());
577  pcResult->insert(pcResult->end(), injectedHits_.begin(), injectedHits_.end());
578  e.put(std::move(pcResult), "HcalHits");
579  }
580 
581 #ifdef EDM_ML_DEBUG
582  edm::LogVerbatim("HcalSim") << "\n========> HcalDigitizer e.put\n";
583 #endif
584 }
585 
587  checkGeometry(es);
588 
589  if (agingFlagHB) {
591  }
592  if (agingFlagHE) {
594  }
595 
597 
600  theHOAmplifier->setTimeSlew(hcalTimeSlew_delay_);
601  theZDCAmplifier->setTimeSlew(hcalTimeSlew_delay_);
602 }
603 
607 
608  if (theHBHEResponse)
609  theHBHEResponse->setGeometry(theGeometry);
611  theHBHESiPMResponse->setGeometry(theGeometry);
612  if (theHOResponse)
613  theHOResponse->setGeometry(theGeometry);
614  if (theHOSiPMResponse)
615  theHOSiPMResponse->setGeometry(theGeometry);
616  theHFResponse->setGeometry(theGeometry);
617  theHFQIE10Response->setGeometry(theGeometry);
618  theZDCResponse->setGeometry(theGeometry);
619  if (theRelabeller)
620  theRelabeller->setGeometry(theRecNumber);
621 
622  // See if it's been updated
623  bool check1 = theGeometryWatcher_.check(eventSetup);
624  bool check2 = theRecNumberWatcher_.check(eventSetup);
625  if (check1 or check2) {
627  }
628 }
629 
631  const std::vector<DetId> &hbCells = theGeometry->getValidDetIds(DetId::Hcal, HcalBarrel);
632  const std::vector<DetId> &heCells = theGeometry->getValidDetIds(DetId::Hcal, HcalEndcap);
633  const std::vector<DetId> &hoCells = theGeometry->getValidDetIds(DetId::Hcal, HcalOuter);
634  const std::vector<DetId> &hfCells = theGeometry->getValidDetIds(DetId::Hcal, HcalForward);
635  const std::vector<DetId> &zdcCells = theGeometry->getValidDetIds(DetId::Calo, HcalZDCDetId::SubdetectorId);
636  // const std::vector<DetId>& hcalTrigCells =
637  // geometry->getValidDetIds(DetId::Hcal, HcalTriggerTower); const
638  // std::vector<DetId>& hcalCalib = geometry->getValidDetIds(DetId::Calo,
639  // HcalCastorDetId::SubdetectorId);
640  // edm::LogVerbatim("HcalSim") <<"HcalDigitizer::CheckGeometry number of cells: << zdcCells.size();
641  if (zdcCells.empty())
642  zdcgeo = false;
643  if (hbCells.empty() && heCells.empty())
644  hbhegeo = false;
645  if (hoCells.empty())
646  hogeo = false;
647  if (hfCells.empty())
648  hfgeo = false;
649  // combine HB & HE
650 
651  hbheCells = hbCells;
652  if (!killHE_) {
653  hbheCells.insert(hbheCells.end(), heCells.begin(), heCells.end());
654  }
655  // handle mixed QIE8/11 scenario in HBHE
659 
660  if (theHOSiPMDigitizer) {
661  buildHOSiPMCells(hoCells, eventSetup);
662  if (theHOSiPMResponse)
663  theHOSiPMResponse->setDetIds(hoCells);
664  }
665 
666  // handle mixed QIE8/10 scenario in HF
667  buildHFQIECells(hfCells, eventSetup);
668 
669  if (theZDCDigitizer)
670  theZDCDigitizer->setDetIds(zdcCells);
671 
672  // fill test hits collection if desired and empty
673  if (injectTestHits_ && injectedHits_.empty() && !injectedHitsCells_.empty() && !injectedHitsEnergy_.empty()) {
674  // make list of specified cells if desired
675  std::vector<DetId> testCells;
676  if (injectedHitsCells_.size() >= 4) {
677  testCells.reserve(injectedHitsCells_.size() / 4);
678  for (unsigned ic = 0; ic < injectedHitsCells_.size(); ic += 4) {
679  if (ic + 4 > injectedHitsCells_.size())
680  break;
681  testCells.push_back(HcalDetId((HcalSubdetector)injectedHitsCells_[ic],
682  injectedHitsCells_[ic + 1],
683  injectedHitsCells_[ic + 2],
684  injectedHitsCells_[ic + 3]));
685  }
686  } else {
687  int testSubdet = injectedHitsCells_[0];
688  if (testSubdet == HcalBarrel)
689  testCells = hbCells;
690  else if (testSubdet == HcalEndcap)
691  testCells = heCells;
692  else if (testSubdet == HcalForward)
693  testCells = hfCells;
694  else if (testSubdet == HcalOuter)
695  testCells = hoCells;
696  else
697  throw cms::Exception("Configuration") << "Unknown subdet " << testSubdet << " for HCAL test hit injection";
698  }
699  bool useHitTimes = (injectedHitsTime_.size() == injectedHitsEnergy_.size());
700  injectedHits_.reserve(testCells.size() * injectedHitsEnergy_.size());
701  for (unsigned ih = 0; ih < injectedHitsEnergy_.size(); ++ih) {
702  double tmp = useHitTimes ? injectedHitsTime_[ih] : 0.;
703  for (auto &aCell : testCells) {
704  injectedHits_.emplace_back(aCell, injectedHitsEnergy_[ih], tmp);
705  }
706  }
707  }
708 }
709 
710 void HcalDigitizer::buildHFQIECells(const std::vector<DetId> &allCells, const edm::EventSetup &eventSetup) {
711  // if results are already cached, no need to look again
712  if (!theHFQIE8DetIds.empty() || !theHFQIE10DetIds.empty())
713  return;
714 
715  // get the QIETypes
716  // intentional copy
717  HcalQIETypes qieTypes = eventSetup.getData(qieTypesToken_);
718  if (qieTypes.topo() == nullptr) {
719  qieTypes.setTopo(&eventSetup.getData(topoToken_));
720  }
721 
722  for (std::vector<DetId>::const_iterator detItr = allCells.begin(); detItr != allCells.end(); ++detItr) {
723  HcalQIENum qieType = HcalQIENum(qieTypes.getValues(*detItr)->getValue());
724  if (qieType == QIE8) {
725  theHFQIE8DetIds.push_back(*detItr);
726  } else if (qieType == QIE10) {
727  theHFQIE10DetIds.push_back(*detItr);
728  } else { // default is QIE8
729  theHFQIE8DetIds.push_back(*detItr);
730  }
731  }
732 
733  if (!theHFQIE8DetIds.empty())
734  theHFDigitizer->setDetIds(theHFQIE8DetIds);
735  else {
736  theHFDigitizer.reset();
737  }
738 
739  if (!theHFQIE10DetIds.empty())
741  else {
742  theHFQIE10Digitizer.reset();
743  }
744 }
745 
746 void HcalDigitizer::buildHBHEQIECells(const std::vector<DetId> &allCells, const edm::EventSetup &eventSetup) {
747  // if results are already cached, no need to look again
748  if (!theHBHEQIE8DetIds.empty() || !theHBHEQIE11DetIds.empty())
749  return;
750 
751  // get the QIETypes
752  // intentional copy
753  HcalQIETypes qieTypes = eventSetup.getData(qieTypesToken_);
754  if (qieTypes.topo() == nullptr) {
755  qieTypes.setTopo(&eventSetup.getData(topoToken_));
756  }
757 
758  for (std::vector<DetId>::const_iterator detItr = allCells.begin(); detItr != allCells.end(); ++detItr) {
759  HcalQIENum qieType = HcalQIENum(qieTypes.getValues(*detItr)->getValue());
760  if (qieType == QIE8) {
761  theHBHEQIE8DetIds.push_back(*detItr);
762  } else if (qieType == QIE11) {
763  theHBHEQIE11DetIds.push_back(*detItr);
764  } else { // default is QIE8
765  theHBHEQIE8DetIds.push_back(*detItr);
766  }
767  }
768 
769  if (!theHBHEQIE8DetIds.empty())
771  else {
772  theHBHEDigitizer.reset();
773  }
774 
775  if (!theHBHEQIE11DetIds.empty())
777  else {
778  theHBHEQIE11Digitizer.reset();
779  }
780 
781  if (!theHBHEQIE8DetIds.empty() && !theHBHEQIE11DetIds.empty()) {
784  }
785 }
786 
787 void HcalDigitizer::buildHOSiPMCells(const std::vector<DetId> &allCells, const edm::EventSetup &eventSetup) {
788  // all HPD
789 
790  if (theHOSiPMCode == 0) {
791  theHODigitizer->setDetIds(allCells);
792  } else if (theHOSiPMCode == 1) {
793  theHOSiPMDigitizer->setDetIds(allCells);
794  // FIXME pick Zecotek or hamamatsu?
795  } else if (theHOSiPMCode == 2) {
796  std::vector<HcalDetId> zecotekDetIds, hamamatsuDetIds;
797 
798  // intentional copy
799  HcalMCParams mcParams = eventSetup.getData(mcParamsToken_);
800  if (mcParams.topo() == nullptr) {
801  mcParams.setTopo(&eventSetup.getData(topoToken_));
802  }
803 
804  for (std::vector<DetId>::const_iterator detItr = allCells.begin(); detItr != allCells.end(); ++detItr) {
805  int shapeType = mcParams.getValues(*detItr)->signalShape();
806  if (shapeType == HcalShapes::ZECOTEK) {
807  zecotekDetIds.emplace_back(*detItr);
808  theHOSiPMDetIds.push_back(*detItr);
809  } else if (shapeType == HcalShapes::HAMAMATSU) {
810  hamamatsuDetIds.emplace_back(*detItr);
811  theHOSiPMDetIds.push_back(*detItr);
812  } else {
813  theHOHPDDetIds.push_back(*detItr);
814  }
815  }
816 
817  if (!theHOHPDDetIds.empty())
818  theHODigitizer->setDetIds(theHOHPDDetIds);
819  else {
820  theHODigitizer.reset();
821  }
822 
823  if (!theHOSiPMDetIds.empty())
825  else {
826  theHOSiPMDigitizer.reset();
827  }
828 
829  if (!theHOHPDDetIds.empty() && !theHOSiPMDetIds.empty()) {
832  }
833 
834  theParameterMap.setHOZecotekDetIds(zecotekDetIds);
835  theParameterMap.setHOHamamatsuDetIds(hamamatsuDetIds);
836 
837  // make sure we don't got through this exercise again
838  theHOSiPMCode = -2;
839  }
840 }
841 
842 void HcalDigitizer::darkening(std::vector<PCaloHit> &hcalHits) {
843  for (unsigned int ii = 0; ii < hcalHits.size(); ++ii) {
844  uint32_t tmpId = hcalHits[ii].id();
845  int det, z, depth, ieta, phi, lay;
846  HcalTestNumbering::unpackHcalIndex(tmpId, det, z, depth, ieta, phi, lay);
847 
848  bool darkened = false;
849  float dweight = 1.;
850 
851  if (det == int(HcalBarrel) && m_HBDarkening) {
852  // HB darkening
853  dweight = m_HBDarkening->degradation(deliveredLumi, ieta, lay);
854  darkened = true;
855  } else if (det == int(HcalEndcap) && m_HEDarkening) {
856  // HE darkening
857  dweight = m_HEDarkening->degradation(deliveredLumi, ieta, lay);
858  darkened = true;
859  } else if (det == int(HcalForward) && m_HFRecalibration) {
860  // HF darkening - approximate: invert recalibration factor
861  dweight = 1.0 / m_HFRecalibration->getCorr(ieta, depth, deliveredLumi);
862  darkened = true;
863  }
864 
865  // reset hit energy
866  if (darkened)
867  hcalHits[ii].setEnergy(hcalHits[ii].energy() * dweight);
868  }
869 }
HFHitFilter theHFQIE10HitFilter
std::vector< DetId > theHFQIE10DetIds
Log< level::Info, true > LogVerbatim
std::vector< int > injectedHitsCells_
HBHEHitFilter theHBHEHitFilter
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
constexpr bool oldFormat() const
Definition: HcalDetId.h:139
std::vector< PCaloHit > PCaloHitContainer
void setQIE10NoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
std::vector< CaloSamples > CaloSamplesCollection
Definition: CaloSamples.h:99
const HcalDDDRecConstants * theRecNumber
void setParameterMap(HcalSimParameterMap *map)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::unique_ptr< CaloHitResponse > theHFResponse
const edm::ESGetToken< HcalDDDRecConstants, HcalRecNumberingRecord > theRecNumberToken
Definition: HcalDigitizer.h:96
std::unique_ptr< HcalTimeSlewSim > theTimeSlewSim
std::unique_ptr< HcalAmplifier > theHBHEAmplifier
int getValue() const
Definition: HcalQIEType.h:19
virtual ~HcalDigitizer()
std::unique_ptr< HBHEDigitizer > theHBHEDigitizer
void buildHFQIECells(const std::vector< DetId > &allCells, const edm::EventSetup &eventSetup)
Definition: HcalQIENum.h:4
float degradation(float intlumi, int ieta, int lay) const
HOHitFilter theHOHitFilter
std::unique_ptr< ZDCDigitizer > theZDCDigitizer
void darkening(std::vector< PCaloHit > &hcalHits)
std::unique_ptr< HcalHitRelabeller > theRelabeller
HFHitFilter theHFHitFilter
void checkGeometry(const edm::EventSetup &eventSetup)
Log< level::Error, false > LogError
void buildHBHEQIECells(const std::vector< DetId > &allCells, const edm::EventSetup &eventSetup)
std::unique_ptr< HFDigitizer > theHFDigitizer
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > theGeometryToken
Definition: HcalDigitizer.h:95
std::unique_ptr< HODigitizer > theHODigitizer
std::unique_ptr< HcalAmplifier > theHFAmplifier
edm::SortedCollection< ZDCDataFrame > ZDCDigiCollection
edm::ESGetToken< HBHEDarkening, HBHEDarkeningRecord > m_HBDarkeningToken
Definition: HcalDigitizer.h:92
const Item * getValues(DetId fId, bool throwOnFail=true) const
std::unique_ptr< HcalCoderFactory > theCoderFactory
HcalDataFrameContainer< QIE10DataFrame > QIE10DigiCollection
void initializeEvent(edm::Event const &e, edm::EventSetup const &c)
std::unique_ptr< HcalElectronicsSim > theZDCElectronicsSim
std::unique_ptr< HFRecalibration > m_HFRecalibration
edm::EDGetTokenT< std::vector< PCaloHit > > zdcToken_
void accumulate(edm::Event const &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *)
constexpr uint32_t newForm() const
Definition: HcalDetId.h:199
std::vector< double > injectedHitsEnergy_
void setQIE11NoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
void setHOZecotekDetIds(const std::vector< HcalDetId > &ids)
std::unique_ptr< HcalElectronicsSim > theHFQIE10ElectronicsSim
const HBHEDarkening * m_HEDarkening
std::unique_ptr< HcalElectronicsSim > theHFElectronicsSim
const HcalTopology * topo() const
edm::ESWatcher< HcalRecNumberingRecord > theRecNumberWatcher_
edm::SortedCollection< HODataFrame > HODigiCollection
std::string hitsProducer_
double deliveredLumi
void setElectronicsSim(HcalElectronicsSim *electronicsSim)
std::unique_ptr< HODigitizer > theHOSiPMDigitizer
Creates electronics signals from hits.
const HcalTimeSlew * hcalTimeSlew_delay_
void setHFNoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
std::unique_ptr< HcalAmplifier > theHBHEQIE11Amplifier
void buildHOSiPMCells(const std::vector< DetId > &allCells, const edm::EventSetup &eventSetup)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
ZDCHitFilter theZDCHitFilter
std::unique_ptr< HcalSiPMHitResponse > theHBHESiPMResponse
HcalSubdetector
Definition: HcalAssistant.h:31
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
std::vector< double > injectedHitsTime_
HcalQIENum
Definition: HcalQIENum.h:4
static const size_type MAXSAMPLES
HcalShapes theShapes
std::unique_ptr< CaloHitResponse > theHBHEResponse
std::vector< DetId > theHOHPDDetIds
edm::EDGetTokenT< std::vector< PCaloHit > > hcalToken_
unsigned int signalShape() const
Definition: HcalMCParam.h:38
edm::ESGetToken< HBHEDarkening, HBHEDarkeningRecord > m_HEDarkeningToken
Definition: HcalDigitizer.h:93
HcalSimParameterMap theParameterMap
std::unique_ptr< HcalElectronicsSim > theHBHEElectronicsSim
ii
Definition: cuy.py:589
virtual bool valid(const HcalZDCDetId &id) const
Definition: ZdcTopology.cc:30
const edm::ESGetToken< HcalQIETypes, HcalQIETypesRcd > qieTypesToken_
Definition: HcalDigitizer.h:97
HcalDataFrameContainer< QIE11DataFrame > QIE11DigiCollection
Log< level::Info, false > LogInfo
const edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > topoToken_
Definition: HcalDigitizer.h:90
std::unique_ptr< CaloHitResponse > theHOResponse
std::unique_ptr< HPDIonFeedbackSim > theIonFeedback
Definition: DetId.h:17
std::unique_ptr< HcalAmplifier > theZDCAmplifier
edm::ESWatcher< CaloGeometryRecord > theGeometryWatcher_
Definition: HcalDigitizer.h:99
void setHOHamamatsuDetIds(const std::vector< HcalDetId > &ids)
void setDbService(const HcalDbService *service)
Definition: HcalShapes.h:24
std::vector< DetId > theHFQIE8DetIds
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
std::unique_ptr< QIE11Digitizer > theHBHEQIE11Digitizer
std::unique_ptr< QIE10Digitizer > theHFQIE10Digitizer
void setDetIds(const std::vector< DetId > &detIds)
Definition: HcalHitFilter.h:17
std::vector< DetId > theHOSiPMDetIds
HOHitFilter theHOSiPMHitFilter
std::unique_ptr< CaloHitResponse > theZDCResponse
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
const CaloGeometry * theGeometry
bool isValid() const
Definition: HandleBase.h:70
std::unique_ptr< HcalAmplifier > theHOAmplifier
std::unique_ptr< HcalSiPMHitResponse > theHOSiPMResponse
bool validHcal(const HcalDetId &id) const
std::vector< DetId > hbheCells
const HBHEDarkening * m_HBDarkening
HLT enums.
void accumulateCaloHits(edm::Handle< std::vector< PCaloHit >> const &hcalHits, edm::Handle< std::vector< PCaloHit >> const &zdcHits, int bunchCrossing, CLHEP::HepRandomEngine *, const HcalTopology *h, const ZdcTopology *z)
std::unique_ptr< HcalElectronicsSim > theHOElectronicsSim
HBHEHitFilter theHBHEQIE11HitFilter
void setZDCNoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
std::vector< DetId > theHBHEQIE8DetIds
std::unique_ptr< CaloHitResponse > theHFQIE10Response
edm::SortedCollection< HFDataFrame > HFDigiCollection
edm::ESGetToken< HcalMCParams, HcalMCParamsRcd > mcParamsToken_
Definition: HcalDigitizer.h:98
const edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
Definition: HcalDigitizer.h:89
static constexpr int32_t SubdetectorId
Definition: HcalZDCDetId.h:35
std::unique_ptr< HcalAmplifier > theHFQIE10Amplifier
std::unique_ptr< HcalElectronicsSim > theHBHEQIE11ElectronicsSim
void updateGeometry(const edm::EventSetup &eventSetup)
tmp
align.sh
Definition: createJobs.py:716
void setHBHENoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
void setup(const edm::EventSetup &es)
const edm::ESGetToken< ZdcTopology, HcalRecNumberingRecord > topoZToken_
Definition: HcalDigitizer.h:91
std::vector< DetId > theHBHEQIE11DetIds
def move(src, dest)
Definition: eostools.py:511
void setTopo(const HcalTopology *topo)
void setDbService(const HcalDbService *service)
edm::SortedCollection< HBHEDataFrame > HBHEDigiCollection
std::vector< PCaloHit > injectedHits_
void finalizeEvent(edm::Event &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *)
void setHONoiseSignalGenerator(HcalBaseSignalGenerator *noiseGenerator)
HcalDigitizer(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
const edm::ESGetToken< HcalTimeSlew, HcalTimeSlewRecord > hcalTimeSlew_delay_token_
Definition: HcalDigitizer.h:94