-
Notifications
You must be signed in to change notification settings - Fork 4.4k
/
Copy pathHcalDigiProducer.cc
100 lines (83 loc) · 3.3 KB
/
HcalDigiProducer.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
#include "SimCalorimetry/HcalSimProducers/interface/HcalDigiProducer.h"
#include "FWCore/Framework/interface/ProducerBase.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Utilities/interface/RandomNumberGenerator.h"
#include "FWCore/Utilities/interface/StreamID.h"
HcalDigiProducer::HcalDigiProducer(edm::ParameterSet const& pset, edm::ProducerBase& mixMod, edm::ConsumesCollector& iC) :
DigiAccumulatorMixMod(),
theDigitizer_(pset, iC) {
mixMod.produces<HBHEDigiCollection>();
mixMod.produces<HODigiCollection>();
mixMod.produces<HFDigiCollection>();
mixMod.produces<ZDCDigiCollection>();
mixMod.produces<QIE10DigiCollection>("HFQIE10DigiCollection");
mixMod.produces<QIE11DigiCollection>("HBHEQIE11DigiCollection");
if(pset.getParameter<bool>("debugCaloSamples")){
mixMod.produces<CaloSamplesCollection>("HcalSamples");
}
if(pset.getParameter<bool>("injectTestHits")){
mixMod.produces<edm::PCaloHitContainer>("HcalHits");
}
}
HcalDigiProducer::HcalDigiProducer(edm::ParameterSet const& pset, edm::ConsumesCollector& iC) :
DigiAccumulatorMixMod(),
theDigitizer_(pset, iC) {
}
void
HcalDigiProducer::initializeEvent(edm::Event const& event, edm::EventSetup const& es) {
theDigitizer_.initializeEvent(event, es);
}
void
HcalDigiProducer::finalizeEvent(edm::Event& event, edm::EventSetup const& es) {
theDigitizer_.finalizeEvent(event, es, randomEngine(event.streamID()));
}
void
HcalDigiProducer::accumulate(edm::Event const& event, edm::EventSetup const& es) {
theDigitizer_.accumulate(event, es, randomEngine(event.streamID()));
}
void
HcalDigiProducer::accumulate(PileUpEventPrincipal const& event, edm::EventSetup const& es, edm::StreamID const& streamID) {
theDigitizer_.accumulate(event, es, randomEngine(streamID));
}
void
HcalDigiProducer::beginRun(edm::Run const&, edm::EventSetup const& es) {}
void
HcalDigiProducer::endRun(edm::Run const&, edm::EventSetup const&) {}
void
HcalDigiProducer::setHBHENoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setHBHENoiseSignalGenerator(noiseGenerator);
}
void
HcalDigiProducer::setHFNoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setHFNoiseSignalGenerator(noiseGenerator);
}
void
HcalDigiProducer::setHONoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setHONoiseSignalGenerator(noiseGenerator);
}
void
HcalDigiProducer::setZDCNoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setZDCNoiseSignalGenerator(noiseGenerator);
}
void
HcalDigiProducer::setQIE10NoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setQIE10NoiseSignalGenerator(noiseGenerator);
}
void
HcalDigiProducer::setQIE11NoiseSignalGenerator(HcalBaseSignalGenerator * noiseGenerator) {
theDigitizer_.setQIE11NoiseSignalGenerator(noiseGenerator);
}
CLHEP::HepRandomEngine* HcalDigiProducer::randomEngine(edm::StreamID const& streamID) {
unsigned int index = streamID.value();
if(index >= randomEngines_.size()) {
randomEngines_.resize(index + 1, nullptr);
}
CLHEP::HepRandomEngine* ptr = randomEngines_[index];
if(!ptr) {
edm::Service<edm::RandomNumberGenerator> rng;
ptr = &rng->getEngine(streamID);
randomEngines_[index] = ptr;
}
return ptr;
}