CMS 3D CMS Logo

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