CMS 3D CMS Logo

DataMixingHcalWorker.cc
Go to the documentation of this file.
1 // File: DataMixingHcalWorker.cc
2 // Description: see DataMixingHcalWorker.h
3 // Author: Mike Hildreth, University of Notre Dame
4 //
5 //--------------------------------------------
6 
9 #include <map>
10 #include <memory>
11 //
12 //
13 #include "DataMixingHcalWorker.h"
14 
15 using namespace std;
16 
17 namespace edm {
18 
19  // Virtual constructor
20 
21  DataMixingHcalWorker::DataMixingHcalWorker() {}
22 
23  // Constructor
24  DataMixingHcalWorker::DataMixingHcalWorker(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
25  : label_(ps.getParameter<std::string>("Label"))
26 
27  {
28  // get the subdetector names
29  // this->getSubdetectorNames(); //something like this may be useful to
30  // check what we are supposed to do...
31 
32  // declare the products to produce
33 
34  // Hcal
35 
36  HBHErechitCollectionSig_ = ps.getParameter<edm::InputTag>("HBHEProducerSig");
37  HOrechitCollectionSig_ = ps.getParameter<edm::InputTag>("HOProducerSig");
38  HFrechitCollectionSig_ = ps.getParameter<edm::InputTag>("HFProducerSig");
39  ZDCrechitCollectionSig_ = ps.getParameter<edm::InputTag>("ZDCrechitCollectionSig");
40 
41  HBHEPileRecHitInputTag_ = ps.getParameter<edm::InputTag>("HBHEPileRecHitInputTag");
42  HOPileRecHitInputTag_ = ps.getParameter<edm::InputTag>("HOPileRecHitInputTag");
43  HFPileRecHitInputTag_ = ps.getParameter<edm::InputTag>("HFPileRecHitInputTag");
44  ZDCPileRecHitInputTag_ = ps.getParameter<edm::InputTag>("ZDCPileRecHitInputTag");
45 
46  HBHERecHitCollectionDM_ = ps.getParameter<std::string>("HBHERecHitCollectionDM");
47  HORecHitCollectionDM_ = ps.getParameter<std::string>("HORecHitCollectionDM");
48  HFRecHitCollectionDM_ = ps.getParameter<std::string>("HFRecHitCollectionDM");
49  ZDCRecHitCollectionDM_ = ps.getParameter<std::string>("ZDCRecHitCollectionDM");
50 
54 
58 
61  }
62 
63  // Virtual destructor needed.
65 
67  // fill in maps of hits
68 
69  LogInfo("DataMixingHcalWorker") << "===============> adding MC signals for " << e.id();
70 
71  // HBHE first
72 
73  Handle<HBHERecHitCollection> pHBHERecHits;
74 
75  const HBHERecHitCollection *HBHERecHits = nullptr;
76 
77  if (e.getByToken(HBHERecHitToken_, pHBHERecHits)) {
78  HBHERecHits = pHBHERecHits.product(); // get a ptr to the product
79  LogDebug("DataMixingHcalWorker") << "total # HBHE rechits: " << HBHERecHits->size();
80  }
81 
82  if (HBHERecHits) {
83  // loop over rechits, storing them in a map so we can add pileup later
84  for (HBHERecHitCollection::const_iterator it = HBHERecHits->begin(); it != HBHERecHits->end(); ++it) {
86 
87 #ifdef DEBUG
88  LogDebug("DataMixingHcalWorker") << "processed HBHERecHit with rawId: " << it->id() << "\n"
89  << " rechit energy: " << it->energy();
90 #endif
91  }
92  }
93 
94  // HO next
95 
96  Handle<HORecHitCollection> pHORecHits;
97 
98  const HORecHitCollection *HORecHits = nullptr;
99 
100  if (e.getByToken(HORecHitToken_, pHORecHits)) {
101  HORecHits = pHORecHits.product(); // get a ptr to the product
102 #ifdef DEBUG
103  LogDebug("DataMixingHcalWorker") << "total # HO rechits: " << HORecHits->size();
104 #endif
105  }
106 
107  if (HORecHits) {
108  // loop over rechits, storing them in a map so we can add pileup later
109  for (HORecHitCollection::const_iterator it = HORecHits->begin(); it != HORecHits->end(); ++it) {
110  HORecHitStorage_.insert(HORecHitMap::value_type((it->id()), *it));
111 
112 #ifdef DEBUG
113  LogDebug("DataMixingHcalWorker") << "processed HORecHit with rawId: " << it->id() << "\n"
114  << " rechit energy: " << it->energy();
115 #endif
116  }
117  }
118 
119  // HF next
120 
121  Handle<HFRecHitCollection> pHFRecHits;
122 
123  const HFRecHitCollection *HFRecHits = nullptr;
124 
125  if (e.getByToken(HFRecHitToken_, pHFRecHits)) {
126  HFRecHits = pHFRecHits.product(); // get a ptr to the product
127 #ifdef DEBUG
128  LogDebug("DataMixingHcalWorker") << "total # HF rechits: " << HFRecHits->size();
129 #endif
130  }
131 
132  if (HFRecHits) {
133  // loop over rechits, storing them in a map so we can add pileup later
134  for (HFRecHitCollection::const_iterator it = HFRecHits->begin(); it != HFRecHits->end(); ++it) {
135  HFRecHitStorage_.insert(HFRecHitMap::value_type((it->id()), *it));
136 
137 #ifdef DEBUG
138  LogDebug("DataMixingHcalWorker") << "processed HFRecHit with rawId: " << it->id() << "\n"
139  << " rechit energy: " << it->energy();
140 #endif
141  }
142  }
143 
144  // ZDC next
145 
146  Handle<ZDCRecHitCollection> pZDCRecHits;
147 
148  const ZDCRecHitCollection *ZDCRecHits = nullptr;
149 
150  if (e.getByToken(ZDCRecHitToken_, pZDCRecHits)) {
151  ZDCRecHits = pZDCRecHits.product(); // get a ptr to the product
152 #ifdef DEBUG
153  LogDebug("DataMixingHcalWorker") << "total # ZDC rechits: " << ZDCRecHits->size();
154 #endif
155  }
156 
157  if (ZDCRecHits) {
158  // loop over rechits, storing them in a map so we can add pileup later
159  for (ZDCRecHitCollection::const_iterator it = ZDCRecHits->begin(); it != ZDCRecHits->end(); ++it) {
160  ZDCRecHitStorage_.insert(ZDCRecHitMap::value_type((it->id()), *it));
161 
162 #ifdef DEBUG
163  LogDebug("DataMixingHcalWorker") << "processed ZDCRecHit with rawId: " << it->id() << "\n"
164  << " rechit energy: " << it->energy();
165 #endif
166  }
167  }
168 
169  } // end of addEMSignals
170 
172  const EventPrincipal *ep,
173  unsigned int eventNr,
174  ModuleCallingContext const *mcc) {
175  LogDebug("DataMixingHcalWorker") << "\n===============> adding pileups from event " << ep->id()
176  << " for bunchcrossing " << bcr;
177 
178  // fill in maps of hits; same code as addSignals, except now applied to the
179  // pileup events
180 
181  // HBHE first
182 
183  std::shared_ptr<Wrapper<HBHERecHitCollection> const> HBHERecHitsPTR =
184  getProductByTag<HBHERecHitCollection>(*ep, HBHEPileRecHitInputTag_, mcc);
185 
186  if (HBHERecHitsPTR) {
187  const HBHERecHitCollection *HBHERecHits = const_cast<HBHERecHitCollection *>(HBHERecHitsPTR->product());
188 
189  LogDebug("DataMixingEMWorker") << "total # HBHE rechits: " << HBHERecHits->size();
190 
191  // loop over rechits, adding these to the existing maps
192  for (HBHERecHitCollection::const_iterator it = HBHERecHits->begin(); it != HBHERecHits->end(); ++it) {
194 
195 #ifdef DEBUG
196  LogDebug("DataMixingEMWorker") << "processed HBHERecHit with rawId: " << it->id().rawId() << "\n"
197  << " rechit energy: " << it->energy();
198 #endif
199  }
200  }
201 
202  // HO Next
203 
204  std::shared_ptr<Wrapper<HORecHitCollection> const> HORecHitsPTR =
205  getProductByTag<HORecHitCollection>(*ep, HOPileRecHitInputTag_, mcc);
206 
207  if (HORecHitsPTR) {
208  const HORecHitCollection *HORecHits = const_cast<HORecHitCollection *>(HORecHitsPTR->product());
209 
210  LogDebug("DataMixingEMWorker") << "total # HO rechits: " << HORecHits->size();
211 
212  // loop over rechits, adding these to the existing maps
213  for (HORecHitCollection::const_iterator it = HORecHits->begin(); it != HORecHits->end(); ++it) {
214  HORecHitStorage_.insert(HORecHitMap::value_type((it->id()), *it));
215 
216 #ifdef DEBUG
217  LogDebug("DataMixingEMWorker") << "processed HORecHit with rawId: " << it->id().rawId() << "\n"
218  << " rechit energy: " << it->energy();
219 #endif
220  }
221  }
222 
223  // HF Next
224 
225  std::shared_ptr<Wrapper<HFRecHitCollection> const> HFRecHitsPTR =
226  getProductByTag<HFRecHitCollection>(*ep, HFPileRecHitInputTag_, mcc);
227 
228  if (HFRecHitsPTR) {
229  const HFRecHitCollection *HFRecHits = const_cast<HFRecHitCollection *>(HFRecHitsPTR->product());
230 
231  LogDebug("DataMixingEMWorker") << "total # HF rechits: " << HFRecHits->size();
232 
233  // loop over rechits, adding these to the existing maps
234  for (HFRecHitCollection::const_iterator it = HFRecHits->begin(); it != HFRecHits->end(); ++it) {
235  HFRecHitStorage_.insert(HFRecHitMap::value_type((it->id()), *it));
236 
237 #ifdef DEBUG
238  LogDebug("DataMixingEMWorker") << "processed HFRecHit with rawId: " << it->id().rawId() << "\n"
239  << " rechit energy: " << it->energy();
240 #endif
241  }
242  }
243 
244  // ZDC Next
245 
246  std::shared_ptr<Wrapper<ZDCRecHitCollection> const> ZDCRecHitsPTR =
247  getProductByTag<ZDCRecHitCollection>(*ep, ZDCPileRecHitInputTag_, mcc);
248 
249  if (ZDCRecHitsPTR) {
250  const ZDCRecHitCollection *ZDCRecHits = const_cast<ZDCRecHitCollection *>(ZDCRecHitsPTR->product());
251 
252  LogDebug("DataMixingEMWorker") << "total # ZDC rechits: " << ZDCRecHits->size();
253 
254  // loop over rechits, adding these to the existing maps
255  for (ZDCRecHitCollection::const_iterator it = ZDCRecHits->begin(); it != ZDCRecHits->end(); ++it) {
256  ZDCRecHitStorage_.insert(ZDCRecHitMap::value_type((it->id()), *it));
257 
258 #ifdef DEBUG
259  LogDebug("DataMixingEMWorker") << "processed ZDCRecHit with rawId: " << it->id().rawId() << "\n"
260  << " rechit energy: " << it->energy();
261 #endif
262  }
263  }
264  }
265 
267  // collection of rechits to put in the event
268  std::unique_ptr<HBHERecHitCollection> HBHErechits(new HBHERecHitCollection);
269  std::unique_ptr<HORecHitCollection> HOrechits(new HORecHitCollection);
270  std::unique_ptr<HFRecHitCollection> HFrechits(new HFRecHitCollection);
271  std::unique_ptr<ZDCRecHitCollection> ZDCrechits(new ZDCRecHitCollection);
272 
273  // loop over the maps we have, re-making individual hits or rechits if
274  // necessary.
275  DetId formerID = 0;
276  DetId currentID;
277  float ESum = 0.;
278  float HBTime = 0.;
279 
280  // HB first...
281 
282  HBHERecHitMap::const_iterator iHBchk;
283 
284  for (HBHERecHitMap::const_iterator iHB = HBHERecHitStorage_.begin(); iHB != HBHERecHitStorage_.end(); ++iHB) {
285  currentID = iHB->first;
286 
287  if (currentID == formerID) { // we have to add these rechits together
288 
289  ESum += (iHB->second).energy();
290 
291  } else {
292  if (formerID > 0) {
293  // cutoff for ESum?
294  HBHERecHit aHit(formerID, ESum, HBTime);
295  HBHErechits->push_back(aHit);
296  }
297  // save pointers for next iteration
298  formerID = currentID;
299  ESum = (iHB->second).energy();
300  HBTime = (iHB->second).time(); // take time of first hit in sequence - is this ok?
301  }
302 
303  iHBchk = iHB;
304  if ((++iHBchk) == HBHERecHitStorage_.end()) { // make sure not to lose the last one
305  HBHERecHit aHit(formerID, ESum, HBTime);
306  HBHErechits->push_back(aHit);
307  }
308  }
309 
310  // HO next...
311 
312  // loop over the maps we have, re-making individual hits or rechits if
313  // necessary.
314  formerID = 0;
315  ESum = 0.;
316  float HOTime = 0.;
317 
318  HORecHitMap::const_iterator iHOchk;
319 
320  for (HORecHitMap::const_iterator iHO = HORecHitStorage_.begin(); iHO != HORecHitStorage_.end(); ++iHO) {
321  currentID = iHO->first;
322 
323  if (currentID == formerID) { // we have to add these rechits together
324 
325  ESum += (iHO->second).energy();
326 
327  } else {
328  if (formerID > 0) {
329  // cutoff for ESum?
330  HORecHit aHit(formerID, ESum, HOTime);
331  HOrechits->push_back(aHit);
332  }
333  // save pointers for next iteration
334  formerID = currentID;
335  ESum = (iHO->second).energy();
336  HOTime = (iHO->second).time(); // take time of first hit in sequence - is this ok?
337  }
338 
339  iHOchk = iHO;
340  if ((++iHOchk) == HORecHitStorage_.end()) { // make sure not to lose the last one
341  HORecHit aHit(formerID, ESum, HOTime);
342  HOrechits->push_back(aHit);
343  }
344  }
345 
346  // HF next...
347 
348  // loop over the maps we have, re-making individual hits or rechits if
349  // necessary.
350  formerID = 0;
351  ESum = 0.;
352  float HFTime = 0.;
353 
354  HFRecHitMap::const_iterator iHFchk;
355 
356  for (HFRecHitMap::const_iterator iHF = HFRecHitStorage_.begin(); iHF != HFRecHitStorage_.end(); ++iHF) {
357  currentID = iHF->first;
358 
359  if (currentID == formerID) { // we have to add these rechits together
360 
361  ESum += (iHF->second).energy();
362 
363  } else {
364  if (formerID > 0) {
365  // cutoff for ESum?
366  HFRecHit aHit(formerID, ESum, HFTime);
367  HFrechits->push_back(aHit);
368  }
369  // save pointers for next iteration
370  formerID = currentID;
371  ESum = (iHF->second).energy();
372  HFTime = (iHF->second).time(); // take time of first hit in sequence - is this ok?
373  }
374 
375  iHFchk = iHF;
376  if ((++iHFchk) == HFRecHitStorage_.end()) { // make sure not to lose the last one
377  HFRecHit aHit(formerID, ESum, HBTime);
378  HFrechits->push_back(aHit);
379  }
380  }
381 
382  // ZDC next...
383 
384  // loop over the maps we have, re-making individual hits or rechits if
385  // necessary.
386  formerID = 0;
387  ESum = 0.;
388  float ZDCTime = 0.;
389  float lowGainEnergy = 0;
390  ZDCRecHit ZOldHit;
391 
392  ZDCRecHitMap::const_iterator iZDCchk;
393 
394  for (ZDCRecHitMap::const_iterator iZDC = ZDCRecHitStorage_.begin(); iZDC != ZDCRecHitStorage_.end(); ++iZDC) {
395  currentID = iZDC->first;
396 
397  if (currentID == formerID) { // we have to add these rechits together
398 
399  ESum += (iZDC->second).energy();
400 
401  } else {
402  if (formerID > 0) {
403  // cutoff for ESum?
404  ZDCRecHit aHit(formerID, ESum, ZDCTime, lowGainEnergy);
405  ZDCrechits->push_back(aHit);
406  }
407  // save pointers for next iteration
408  formerID = currentID;
409  ESum = (iZDC->second).energy();
410  lowGainEnergy = (iZDC->second).lowGainEnergy();
411  ZDCTime = (iZDC->second).time(); // take time of first hit in sequence - is this ok?
412  }
413 
414  iZDCchk = iZDC;
415  if ((++iZDCchk) == ZDCRecHitStorage_.end()) { // make sure not to lose the last one
416  ZDCRecHit aHit(formerID, ESum, HBTime, lowGainEnergy);
417  ZDCrechits->push_back(aHit);
418  }
419  }
420 
421  // done merging
422 
423  // put the collection of recunstructed hits in the event
424  LogInfo("DataMixingHcalWorker") << "total # HBHE Merged rechits: " << HBHErechits->size();
425  LogInfo("DataMixingHcalWorker") << "total # HO Merged rechits: " << HOrechits->size();
426  LogInfo("DataMixingHcalWorker") << "total # HF Merged rechits: " << HFrechits->size();
427  LogInfo("DataMixingHcalWorker") << "total # ZDC Merged rechits: " << ZDCrechits->size();
428 
429  e.put(std::move(HBHErechits), HBHERecHitCollectionDM_);
430  e.put(std::move(HOrechits), HORecHitCollectionDM_);
431  e.put(std::move(HFrechits), HFRecHitCollectionDM_);
432  e.put(std::move(ZDCrechits), ZDCRecHitCollectionDM_);
433 
434  // clear local storage after this event
435  HBHERecHitStorage_.clear();
436  HORecHitStorage_.clear();
437  HFRecHitStorage_.clear();
438  ZDCRecHitStorage_.clear();
439  }
440 
441 } // namespace edm
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void addHcalSignals(const edm::Event &e)
void addHcalPileups(const int bcr, const edm::EventPrincipal *, unsigned int EventId, ModuleCallingContext const *)
size_type size() const
T const * product() const
Definition: Handle.h:70
std::vector< T >::const_iterator const_iterator
edm::EDGetTokenT< ZDCRecHitCollection > ZDCRecHitToken_
edm::EDGetTokenT< ZDCRecHitCollection > ZDCRecHitPToken_
edm::EDGetTokenT< HORecHitCollection > HORecHitToken_
edm::EDGetTokenT< HBHERecHitCollection > HBHERecHitPToken_
edm::EDGetTokenT< HFRecHitCollection > HFRecHitToken_
const_iterator begin() const
const_iterator end() const
Log< level::Info, false > LogInfo
Definition: DetId.h:17
edm::EDGetTokenT< HBHERecHitCollection > HBHERecHitToken_
edm::EDGetTokenT< HFRecHitCollection > HFRecHitPToken_
HLT enums.
edm::EDGetTokenT< HORecHitCollection > HORecHitPToken_
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)