CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalTBDigiProducer.cc
Go to the documentation of this file.
8 
16 
18 
20  edm::ProducesCollector producesCollector,
22  : theParameterMap(new HcalTBSimParameterMap(ps)),
23  paraMap(new HcalSimParameterMap(ps)),
24  theHcalShape(new HcalShape()),
25  theHcalIntegratedShape(new CaloShapeIntegrator(theHcalShape)),
26  theHBHEResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
27  theHOResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
28  theAmplifier(nullptr),
29  theCoderFactory(nullptr),
30  theElectronicsSim(nullptr),
31  theTimeSlewSim(nullptr),
32  theHBHEDigitizer(nullptr),
33  theHODigitizer(nullptr),
34  conditionsToken_(iC.esConsumes()),
35  hcalTimeSlew_delay_token_(iC.esConsumes(edm::ESInputTag("", "HBHE"))),
36  theHBHEHits(),
37  theHOHits(),
38  thisPhaseShift(0) {
39  std::string const instance("simHcalDigis");
40  producesCollector.produces<HBHEDigiCollection>(instance);
41  producesCollector.produces<HODigiCollection>(instance);
42  iC.consumes<std::vector<PCaloHit>>(edm::InputTag("g4SimHits", "HcalHits"));
43 
44  DetId detId(DetId::Hcal, 1);
45  bool syncPhase = (theParameterMap->simParameters(detId)).syncPhase();
46  doPhaseShift = !syncPhase;
47 
50 
51  bool doNoise = ps.getParameter<bool>("doNoise");
52  bool dummy1 = false;
53  bool dummy2 = false; // extra arguments for premixing
54  theAmplifier = new HcalAmplifier(theParameterMap, doNoise, dummy1, dummy2);
57 
58  double minFCToDelay = ps.getParameter<double>("minFCToDelay");
59  bool doTimeSlew = ps.getParameter<bool>("doTimeSlew");
60 
61  hcalTimeSlew_delay_ = nullptr;
62  if (doTimeSlew) {
63  // no time slewing for HF
64  theTimeSlewSim = new HcalTimeSlewSim(theParameterMap, minFCToDelay);
66  }
67 
70 
71  tunePhaseShift = ps.getUntrackedParameter<double>("tunePhaseShiftTB", 1.);
72  ecalTBInfoLabel = ps.getUntrackedParameter<std::string>("EcalTBInfoLabel", "SimEcalTBG4Object");
73  edm::LogInfo("HcalSim") << "HcalTBDigiProducer initialized with doNoise = " << doNoise
74  << ", doTimeSlew = " << doTimeSlew << " and doPhaseShift = " << doPhaseShift
75  << " tunePhasShift = " << tunePhaseShift;
76 
77  if (doPhaseShift) {
79  }
80 }
81 
83  if (theHBHEDigitizer)
84  delete theHBHEDigitizer;
85  if (theHODigitizer)
86  delete theHODigitizer;
87  if (theParameterMap)
88  delete theParameterMap;
89  if (theHcalShape)
90  delete theHcalShape;
93  if (theHBHEResponse)
94  delete theHBHEResponse;
95  if (theHOResponse)
96  delete theHOResponse;
98  delete theElectronicsSim;
99  if (theAmplifier)
100  delete theAmplifier;
101  if (theCoderFactory)
102  delete theCoderFactory;
103  if (theTimeSlewSim)
104  delete theTimeSlewSim;
105 }
106 
108  // get the appropriate gains, noises, & widths for this event
109  const HcalDbService *conditions = &eventSetup.getData(conditionsToken_);
110  theAmplifier->setDbService(conditions);
111  theCoderFactory->setDbService(conditions);
112 
113  // get the correct geometry
114  checkGeometry(eventSetup);
115 
116  // Cache random number engine
118  randomEngine_ = &rng->getEngine(e.streamID());
119 
120  theHBHEHits.clear();
121  theHOHits.clear();
122  if (doPhaseShift) {
123  edm::Handle<PEcalTBInfo> theEcalTBInfo;
124  e.getByLabel(ecalTBInfoLabel, theEcalTBInfo);
125  thisPhaseShift = theEcalTBInfo->phaseShift();
126 
127  DetId detIdHB(DetId::Hcal, 1);
128  setPhaseShift(detIdHB);
129  DetId detIdHO(DetId::Hcal, 3);
130  setPhaseShift(detIdHO);
131  }
132 
134 
136 
139 }
140 
141 void HcalTBDigiProducer::accumulateCaloHits(edm::Handle<std::vector<PCaloHit>> const &hcalHandle, int bunchCrossing) {
142  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate trying to get SimHit";
143 
144  if (hcalHandle.isValid()) {
145  std::vector<PCaloHit> hits = *hcalHandle.product();
146  LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate Hits corrected";
147  theHBHEDigitizer->add(hits, bunchCrossing, randomEngine_);
148  theHODigitizer->add(hits, bunchCrossing, randomEngine_);
149  }
150 }
151 
153  // Step A: Get Inputs, and accumulate digis
154 
155  edm::InputTag hcalTag("g4SimHits", "HcalHits");
157  e.getByLabel(hcalTag, hcalHandle);
158 
159  accumulateCaloHits(hcalHandle, 0);
160 }
161 
163  edm::EventSetup const &,
164  edm::StreamID const &streamID) {
165  // Step A: Get Inputs, and accumulate digis
166 
167  edm::InputTag hcalTag("g4SimHits", "HcalHits");
169  e.getByLabel(hcalTag, hcalHandle);
170 
171  accumulateCaloHits(hcalHandle, e.bunchCrossing());
172 }
173 
175  // Step B: Create empty output
176  std::unique_ptr<HBHEDigiCollection> hbheResult(new HBHEDigiCollection());
177  std::unique_ptr<HODigiCollection> hoResult(new HODigiCollection());
178  LogDebug("HcalSim") << "HcalTBDigiProducer::produce Empty collection created";
179  // Step C: Invoke the algorithm, getting back outputs.
180  theHBHEDigitizer->run(*hbheResult, randomEngine_);
181  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HBHE digis : " << hbheResult->size();
182  theHODigitizer->run(*hoResult, randomEngine_);
183  edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HO digis : " << hoResult->size();
184 
185  // Step D: Put outputs into event
186  std::string const instance("simHcalDigis");
187  e.put(std::move(hbheResult), instance);
188  e.put(std::move(hoResult), instance);
189 
190  randomEngine_ = nullptr; // to prevent access outside event
191 }
192 
194  for (edm::PCaloHitContainer::const_iterator hitItr = hits.begin(); hitItr != hits.end(); ++hitItr) {
195  HcalSubdetector subdet = HcalDetId(hitItr->id()).subdet();
196  if (subdet == HcalBarrel || subdet == HcalEndcap) {
197  theHBHEHits.push_back(*hitItr);
198  } else if (subdet == HcalOuter) {
199  theHOHits.push_back(*hitItr);
200  } else {
201  edm::LogError("HcalSim") << "Bad HcalHit subdetector " << subdet;
202  }
203  }
204 }
205 
207  // see if we need to update
208  if (geometryWatcher_.check(eventSetup)) {
209  theGeometry = &eventSetup.getData(geometryToken_);
210  updateGeometry();
211  }
212 }
213 
217 
218  // Get cells for HB and HE
219  hbheCells.clear();
221  std::vector<DetId> heCells = theGeometry->getValidDetIds(DetId::Hcal, HcalEndcap);
222  // combine HB & HE
223  hbheCells.insert(hbheCells.end(), heCells.begin(), heCells.end());
224 
225  // Get cells for HO
226  hoCells.clear();
228 
229  edm::LogInfo("HcalSim") << "HcalTBDigiProducer update Geometry with " << hbheCells.size() << " cells in HB/HE and "
230  << hoCells.size() << " cells in HO";
231 
233  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HB/HE";
235  LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HO";
236 }
237 
240  if (!parameters.syncPhase()) {
241  int myDet = detId.subdetId();
242  double passPhaseShift = thisPhaseShift + tunePhaseShift;
243  if (myDet <= 2) {
244  theHBHEResponse->setPhaseShift(passPhaseShift);
245  } else {
246  theHOResponse->setPhaseShift(passPhaseShift);
247  }
248  }
249 }
void setPhaseShift(const DetId &detId)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void setGeometry(const CaloGeometry *geometry)
geometry needed for time-of-flight
void setDbService(const HcalDbService *service)
the Producer will probably update this every event
T getUntrackedParameter(std::string const &, T const &) const
void finalizeEvent(edm::Event &e, edm::EventSetup const &c) override
void sortHits(const edm::PCaloHitContainer &hits)
fills the vectors for each subdetector
std::vector< PCaloHit > PCaloHitContainer
std::vector< DetId > hoCells
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
CaloHitResponse * theHBHEResponse
void add(const std::vector< PCaloHit > &hits, int bunchCrossing, CLHEP::HepRandomEngine *engine)
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
static PFTauRenderPlugin instance
HBHEHitFilter theHBHEHitFilter
const CaloSimParameters & simParameters(const DetId &id) const override
edm::ESGetToken< HcalTimeSlew, HcalTimeSlewRecord > hcalTimeSlew_delay_token_
HcalSimParameterMap * paraMap
HODigitizer * theHODigitizer
HcalElectronicsSim * theElectronicsSim
HOHitFilter theHOHitFilter
Log< level::Error, false > LogError
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
HcalTBSimParameterMap * theParameterMap
Main class for Parameters in different subdetectors.
const CaloGeometry * theGeometry
void setPhaseShift(const double &thePhaseShift)
setting the phase shift for asynchronous trigger (e.g. test beams)
CaloTDigitizer< HBHEDigitizerTraits > HBHEDigitizer
std::vector< PCaloHit > theHBHEHits
CaloTDigitizer< HODigitizerTraits > HODigitizer
bool getData(T &iHolder) const
Definition: EventSetup.h:128
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
edm::SortedCollection< HODataFrame > HODigiCollection
CaloVShape * theHcalIntegratedShape
void setHitFilter(const CaloVHitFilter *filter)
Creates electronics signals from hits.
void initializeHits()
HBHEDigitizer * theHBHEDigitizer
void setTimeSlew(const HcalTimeSlew *timeSlew)
Definition: HcalAmplifier.h:38
def move
Definition: eostools.py:511
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
HcalSubdetector
Definition: HcalAssistant.h:31
void checkGeometry(const edm::EventSetup &eventSetup)
shaper for Hcal (not for HF)
Definition: HcalShape.h:15
HcalTBDigiProducer(const edm::ParameterSet &ps, edm::ProducesCollector, edm::ConsumesCollector &iC)
std::vector< DetId > hbheCells
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
void accumulate(edm::Event const &e, edm::EventSetup const &c) override
std::vector< PCaloHit > theHOHits
edm::ESWatcher< CaloGeometryRecord > geometryWatcher_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometryToken_
Log< level::Info, false > LogInfo
Definition: DetId.h:17
CaloHitResponse * theHOResponse
void setTimeSlewSim(HcalTimeSlewSim *timeSlewSim)
Definition: HcalAmplifier.h:35
~HcalTBDigiProducer() override
void setDbService(const HcalDbService *service)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void run(MixCollection< PCaloHit > &, DigiCollection &)
turns hits into digis
HcalAmplifier * theAmplifier
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
std::string ecalTBInfoLabel
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
CLHEP::HepRandomEngine * randomEngine_
HcalCoderFactory * theCoderFactory
CaloVShape * theHcalShape
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
StreamID streamID() const
Definition: Event.h:98
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const HcalTimeSlew * hcalTimeSlew_delay_
void accumulateCaloHits(edm::Handle< std::vector< PCaloHit >> const &hits, int bunchCrossing)
void setDetIds(const std::vector< DetId > &detIds)
edm::SortedCollection< HBHEDataFrame > HBHEDigiCollection
bool syncPhase() const
HcalTimeSlewSim * theTimeSlewSim
#define LogDebug(id)