CMS 3D CMS Logo

DataMixingHcalDigiWorker.cc
Go to the documentation of this file.
1 // File: DataMixingHcalDigiWorker.cc
2 // Description: see DataMixingHcalDigiWorker.h
3 // Author: Mike Hildreth, University of Notre Dame
4 //
5 //--------------------------------------------
6 
7 #include <map>
8 #include <memory>
18 // calibration headers, for future reference
23 //
24 //
26 
27 using namespace std;
28 
29 namespace {
30 
31  typedef std::multimap<DetId, CaloSamples> HcalDigiMap;
32 
33  template <class DIGI>
34  void convertFc2adc (const CaloSamples& fC, const HcalDbService& conditions, DIGI& digi, int capIdOffset = 0) {
35  HcalDetId id = fC.id();
36  const HcalQIECoder* channelCoder = conditions.getHcalCoder (id);
37  const HcalQIEShape* shape = conditions.getHcalShape (channelCoder);
38  HcalCoderDb coder (*channelCoder, *shape);
39  coder.fC2adc(fC, digi, capIdOffset);
40  }
41 
42  template <class DIGI>
43  void convertAdc2fChelper (const DIGI& digi, const HcalDbService& conditions, CaloSamples& fC, const HcalDetId& id) {
44  const HcalCalibrations& calib = conditions.getHcalCalibrations(id);
45  for (int i = 0; i < digi.size(); ++i) {
46  int capId (digi.sample(i).capid());
47  fC[i] -= calib.pedestal (capId);
48  }
49  }
50 
51  template <>
52  void convertAdc2fChelper<QIE10DataFrame> (const QIE10DataFrame& digi, const HcalDbService& conditions, CaloSamples& fC, const HcalDetId& id) {
53  const HcalCalibrations& calib = conditions.getHcalCalibrations(id);
54  for (int i = 0; i < digi.samples(); ++i) {
55  int capId (digi[i].capid());
56  fC[i] -= calib.pedestal (capId);
57  }
58  }
59 
60  template <>
61  void convertAdc2fChelper<QIE11DataFrame> (const QIE11DataFrame& digi, const HcalDbService& conditions, CaloSamples& fC, const HcalDetId& id) {
62  const HcalCalibrations& calib = conditions.getHcalCalibrations(id);
63  for (int i = 0; i < digi.samples(); ++i) {
64  int capId (digi[i].capid());
65  fC[i] -= calib.pedestal (capId);
66  }
67  }
68 
69  template <class DIGI>
70  void convertAdc2fC (const DIGI& digi, const HcalDbService& conditions, bool keepPedestals, CaloSamples& fC) {
71  HcalDetId id (digi.id());
72  const HcalQIECoder* channelCoder = conditions.getHcalCoder (id);
73  const HcalQIEShape* shape = conditions.getHcalShape (channelCoder);
74  HcalCoderDb coder (*channelCoder, *shape);
75  coder.adc2fC(digi, fC);
76  if (!keepPedestals) { // subtract pedestals
77  convertAdc2fChelper(digi,conditions,fC,id);
78  }
79  }
80 
81  template <class DIGIS>
82  void convertHcalDigis (const DIGIS& digis, const HcalDbService& conditions, bool keepPedestals, HcalDigiMap& map) {
83  for (auto digi = digis.begin(); digi != digis.end(); ++digi) {
84  CaloSamples fC;
85  convertAdc2fC (*digi, conditions, keepPedestals, fC);
86  if (!keepPedestals && map.find(digi->id()) == map.end()) {
87  edm::LogWarning("DataMixingHcalDigiWorker")<<"No signal hits found for HCAL cell "<<digi->id()<<" Pedestals may be lost for mixed hit";
88  }
89  map.insert(HcalDigiMap::value_type (digi->id(), fC));
90  }
91  }
92 
93  template <>
94  void convertHcalDigis<QIE10DigiCollection> (const QIE10DigiCollection& digis, const HcalDbService& conditions, bool keepPedestals, HcalDigiMap& map) {
95  for (auto digiItr = digis.begin(); digiItr != digis.end(); ++digiItr) {
96  QIE10DataFrame digi(*digiItr);
97  CaloSamples fC;
98  convertAdc2fC (digi, conditions, keepPedestals, fC);
99  if (!keepPedestals && map.find(digi.id()) == map.end()) {
100  edm::LogWarning("DataMixingHcalDigiWorker")<<"No signal hits found for HCAL cell "<<digi.id()<<" Pedestals may be lost for mixed hit";
101  }
102  map.insert(HcalDigiMap::value_type (digi.id(), fC));
103  }
104  }
105 
106  template <>
107  void convertHcalDigis<QIE11DigiCollection> (const QIE11DigiCollection& digis, const HcalDbService& conditions, bool keepPedestals, HcalDigiMap& map) {
108  for (auto digiItr = digis.begin(); digiItr != digis.end(); ++digiItr) {
109  QIE11DataFrame digi(*digiItr);
110  CaloSamples fC;
111  convertAdc2fC (digi, conditions, keepPedestals, fC);
112  if (!keepPedestals && map.find(digi.id()) == map.end()) {
113  edm::LogWarning("DataMixingHcalDigiWorker")<<"No signal hits found for HCAL cell "<<digi.id()<<" Pedestals may be lost for mixed hit";
114  }
115  map.insert(HcalDigiMap::value_type (digi.id(), fC));
116  }
117  }
118 
119  template <class DIGIS>
120  bool convertSignalHcalDigis (const edm::Event &e, const edm::EDGetTokenT<DIGIS>& token,
121  const HcalDbService& conditions, HcalDigiMap& map) {
122  edm::Handle<DIGIS> digis;
123  if (!e.getByToken (token, digis)) return false;
124  convertHcalDigis (*digis, conditions, true, map); // keep pedestals
125  return true;
126  }
127 
128  template <class DIGIS>
129  bool convertPileupHcalDigis (const edm::EventPrincipal& ep,
130  const edm::InputTag& tag, const edm::ModuleCallingContext* mcc,
131  const HcalDbService& conditions, HcalDigiMap& map) {
132  auto digis = edm::getProductByTag<DIGIS>(ep, tag, mcc);
133  if (!digis) return false;
134  convertHcalDigis (*(digis->product()), conditions, false, map); // subtract pedestals
135  return true;
136  }
137 
138  template <class DIGIS>
139  void buildHcalDigisHelper (DIGIS& digis, const DetId& formerID, CaloSamples& resultSample, const HcalDbService& conditions){
140  // make new digi
141  digis.push_back(typename DIGIS::value_type(formerID));
142  convertFc2adc (resultSample, conditions, digis.back(), 0); // FR guess: simulation starts with 0
143  }
144 
145  template <>
146  void buildHcalDigisHelper<QIE10DigiCollection> (QIE10DigiCollection& digis, const DetId& formerID, CaloSamples& resultSample, const HcalDbService& conditions){
147  // make new digi
148  digis.push_back(formerID.rawId());
149  QIE10DataFrame digi( digis.back() );
150  convertFc2adc (resultSample, conditions, digi, 0); // FR guess: simulation starts with 0
151  }
152 
153  template <>
154  void buildHcalDigisHelper<QIE11DigiCollection> (QIE11DigiCollection& digis, const DetId& formerID, CaloSamples& resultSample, const HcalDbService& conditions){
155  // make new digi
156  digis.push_back(formerID.rawId());
157  QIE11DataFrame digi( digis.back() );
158  convertFc2adc (resultSample, conditions, digi, 0); // FR guess: simulation starts with 0
159  }
160 
161  template <class DIGIS>
162  std::unique_ptr<DIGIS> buildHcalDigis (const HcalDigiMap& map, const HcalDbService& conditions) {
163  std::unique_ptr<DIGIS> digis( new DIGIS );
164  // loop over the maps we have, re-making individual hits or digis if necessary.
165  DetId formerID = 0;
166  CaloSamples resultSample;
167 
168  for(auto hit = map.begin(); hit != map.end(); ++hit) {
169  DetId currentID = hit->first;
170  const CaloSamples& hitSample = hit->second;
171 
172  if (currentID == formerID) { // accumulating hits
173  //loop over digi samples in each CaloSample
174  unsigned int sizenew = (hitSample).size();
175  unsigned int sizeold = resultSample.size();
176  if (sizenew > sizeold) { // extend sample
177  for (unsigned int isamp = sizeold; isamp < sizenew; ++isamp) resultSample[isamp] = 0;
178  resultSample.setSize (sizenew);
179  }
180  for(unsigned int isamp = 0; isamp<sizenew; isamp++) { // add new values
181  resultSample[isamp] += hitSample[isamp]; // for debugging below
182  }
183  }
184  auto hit1 = hit;
185  bool lastEntry = (++hit1 == map.end());
186  if (currentID != formerID || lastEntry) { // store current digi
187  if (formerID>0 || lastEntry) {
188  // make new digi
189  buildHcalDigisHelper(*digis,formerID,resultSample,conditions);
190  }
191  //reset pointers for next iteration
192  formerID = currentID;
193  resultSample = hitSample;
194  }
195  }
196  return digis;
197  }
198 
199 
200 
201 } // namespace {}
202 
203 namespace edm
204 {
205 
206  // Virtual constructor
207 
208  DataMixingHcalDigiWorker::DataMixingHcalDigiWorker() { }
209 
210  // Constructor
211  DataMixingHcalDigiWorker::DataMixingHcalDigiWorker(const edm::ParameterSet& ps, edm::ConsumesCollector && iC) :
212  label_(ps.getParameter<std::string>("Label"))
213  {
214 
215  // get the subdetector names
216  // this->getSubdetectorNames(); //something like this may be useful to check what we are supposed to do...
217 
218  // declare the products to produce
219 
220  // Hcal
221 
222  HBHEdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HBHEdigiCollectionSig");
223  HOdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HOdigiCollectionSig");
224  HFdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HFdigiCollectionSig");
225  ZDCdigiCollectionSig_ = ps.getParameter<edm::InputTag>("ZDCdigiCollectionSig");
226  QIE10digiCollectionSig_ = ps.getParameter<edm::InputTag>("QIE10digiCollectionSig");
227  QIE11digiCollectionSig_ = ps.getParameter<edm::InputTag>("QIE11digiCollectionSig");
228 
229  HBHEPileInputTag_ = ps.getParameter<edm::InputTag>("HBHEPileInputTag");
230  HOPileInputTag_ = ps.getParameter<edm::InputTag>("HOPileInputTag");
231  HFPileInputTag_ = ps.getParameter<edm::InputTag>("HFPileInputTag");
232  ZDCPileInputTag_ = ps.getParameter<edm::InputTag>("ZDCPileInputTag");
233  QIE10PileInputTag_ = ps.getParameter<edm::InputTag>("QIE10PileInputTag");
234  QIE11PileInputTag_ = ps.getParameter<edm::InputTag>("QIE11PileInputTag");
235 
241 
247 
248  DoZDC_ = false;
249  if(ZDCPileInputTag_.label() != "") DoZDC_ = true;
250 
251  if(DoZDC_) {
254  }
255 
256 
257  HBHEDigiCollectionDM_ = ps.getParameter<std::string>("HBHEDigiCollectionDM");
258  HODigiCollectionDM_ = ps.getParameter<std::string>("HODigiCollectionDM");
259  HFDigiCollectionDM_ = ps.getParameter<std::string>("HFDigiCollectionDM");
260  ZDCDigiCollectionDM_ = ps.getParameter<std::string>("ZDCDigiCollectionDM");
261  QIE10DigiCollectionDM_ = ps.getParameter<std::string>("QIE10DigiCollectionDM");
262  QIE11DigiCollectionDM_ = ps.getParameter<std::string>("QIE11DigiCollectionDM");
263 
264 
265  }
266 
267  // Virtual destructor needed.
269  }
270 
272  // Calibration stuff will look like this:
273  // get conditions
274  edm::ESHandle<HcalDbService> conditions;
275  ES.get<HcalDbRecord>().get(conditions);
276 
277 
278  // fill in maps of hits
279 
280  LogInfo("DataMixingHcalDigiWorker")<<"===============> adding MC signals for "<<e.id();
281  convertSignalHcalDigis<HBHEDigiCollection> (e, HBHEDigiToken_, *conditions, HBHEDigiStorage_);
282  convertSignalHcalDigis<HODigiCollection> (e, HODigiToken_, *conditions, HODigiStorage_);
283  convertSignalHcalDigis<HFDigiCollection> (e, HFDigiToken_, *conditions, HFDigiStorage_);
284  convertSignalHcalDigis<QIE10DigiCollection> (e, QIE10DigiToken_, *conditions, QIE10DigiStorage_);
285  convertSignalHcalDigis<QIE11DigiCollection> (e, QIE11DigiToken_, *conditions, QIE11DigiStorage_);
286 
287 
288 
289  // ZDC next
290 
291  if(DoZDC_){
292 
293  Handle< ZDCDigiCollection > pZDCDigis;
294 
295  const ZDCDigiCollection* ZDCDigis = 0;
296 
297  if( e.getByToken( ZDCDigiToken_, pZDCDigis) ) {
298  ZDCDigis = pZDCDigis.product(); // get a ptr to the product
299 #ifdef DEBUG
300  LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
301 #endif
302  }
303 
304 
305  if (ZDCDigis)
306  {
307  // loop over digis, storing them in a map so we can add pileup later
308  for(ZDCDigiCollection::const_iterator it = ZDCDigis->begin(); it != ZDCDigis->end(); ++it) {
309 
310  // calibration, for future reference: (same block for all Hcal types) ZDC is different
311  HcalZDCDetId cell = it->id();
312  // const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
313  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
314  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder); // this one is generic
315  HcalCoderDb coder (*channelCoder, *shape);
316 
317  CaloSamples tool;
318  coder.adc2fC((*it),tool);
319 
320  ZDCDigiStorage_.insert(ZDCDigiMap::value_type( ( it->id() ), tool ));
321 
322 #ifdef DEBUG
323  // Commented out because this does not compile anymore
324  // LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with rawId: "
325  // << it->id() << "\n"
326  // << " digi energy: " << it->energy();
327 #endif
328 
329  }
330  }
331  }
332 
333  } // end of addHCalSignals
334 
335  void DataMixingHcalDigiWorker::addHcalPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr,const edm::EventSetup& ES,
336  ModuleCallingContext const* mcc) {
337 
338  LogDebug("DataMixingHcalDigiWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr;
339 
340  // get conditions
341  edm::ESHandle<HcalDbService> conditions;
342  ES.get<HcalDbRecord>().get(conditions);
343 
344  convertPileupHcalDigis<HBHEDigiCollection> (*ep, HBHEPileInputTag_, mcc, *conditions, HBHEDigiStorage_);
345  convertPileupHcalDigis<HODigiCollection> (*ep, HOPileInputTag_, mcc, *conditions, HODigiStorage_);
346  convertPileupHcalDigis<HFDigiCollection> (*ep, HFPileInputTag_, mcc, *conditions, HFDigiStorage_);
347  convertPileupHcalDigis<QIE10DigiCollection> (*ep, QIE10PileInputTag_, mcc, *conditions, QIE10DigiStorage_);
348  convertPileupHcalDigis<QIE11DigiCollection> (*ep, QIE11PileInputTag_, mcc, *conditions, QIE11DigiStorage_);
349 
350 
351  // ZDC Next
352 
353  if(DoZDC_) {
354 
355 
356  std::shared_ptr<Wrapper<ZDCDigiCollection> const> ZDCDigisPTR = getProductByTag<ZDCDigiCollection>(*ep, ZDCPileInputTag_, mcc);
357 
358  if(ZDCDigisPTR ) {
359 
360  const ZDCDigiCollection* ZDCDigis = const_cast< ZDCDigiCollection * >(ZDCDigisPTR->product());
361 
362  LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
363 
364  // loop over digis, adding these to the existing maps
365  for(ZDCDigiCollection::const_iterator it = ZDCDigis->begin(); it != ZDCDigis->end(); ++it) {
366 
367  // calibration, for future reference: (same block for all Hcal types) ZDC is different
368  HcalZDCDetId cell = it->id();
369  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
370  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder); // this one is generic
371  HcalCoderDb coder (*channelCoder, *shape);
372 
373  CaloSamples tool;
374  coder.adc2fC((*it),tool);
375 
376  ZDCDigiStorage_.insert(ZDCDigiMap::value_type( (it->id()), tool ));
377 
378 #ifdef DEBUG
379  // Commented out because this does not compile anymore
380  // LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with rawId: "
381  // << it->id() << "\n"
382  // << " digi energy: " << it->energy();
383 #endif
384  }
385  }
386  }
387 
388 
389  }
390 
392  edm::ESHandle<HcalDbService> conditions;
393  ES.get<HcalDbRecord>().get(conditions);
394 
395  // collection of digis to put in the event
396  std::unique_ptr< HBHEDigiCollection > HBHEdigis = buildHcalDigis<HBHEDigiCollection> (HBHEDigiStorage_, *conditions);
397  std::unique_ptr< HODigiCollection > HOdigis = buildHcalDigis<HODigiCollection> (HODigiStorage_, *conditions);
398  std::unique_ptr< HFDigiCollection > HFdigis = buildHcalDigis<HFDigiCollection> (HFDigiStorage_, *conditions);
399  std::unique_ptr< QIE10DigiCollection > QIE10digis = buildHcalDigis<QIE10DigiCollection> (QIE10DigiStorage_, *conditions);
400  std::unique_ptr< QIE11DigiCollection > QIE11digis = buildHcalDigis<QIE11DigiCollection> (QIE11DigiStorage_, *conditions);
401  std::unique_ptr< ZDCDigiCollection > ZDCdigis( new ZDCDigiCollection );
402 
403  // loop over the maps we have, re-making individual hits or digis if necessary.
404  DetId formerID = 0;
405  DetId currentID;
406 
407  double fC_new;
408  double fC_old;
409  double fC_sum;
410 
411 
412  // ZDC next...
413 
414  // loop over the maps we have, re-making individual hits or digis if necessary.
415  formerID = 0;
416  CaloSamples ZDC_old;
417 
418  ZDCDigiMap::const_iterator iZDCchk;
419 
420  for(ZDCDigiMap::const_iterator iZDC = ZDCDigiStorage_.begin(); iZDC != ZDCDigiStorage_.end(); ++iZDC) {
421 
422  currentID = iZDC->first;
423 
424  if (currentID == formerID) { // we have to add these digis together
425 
426  //loop over digi samples in each CaloSample
427  unsigned int sizenew = (iZDC->second).size();
428  unsigned int sizeold = ZDC_old.size();
429 
430  unsigned int max_samp = std::max(sizenew, sizeold);
431 
432  CaloSamples ZDC_bigger(currentID,max_samp);
433 
434  bool usenew = false;
435 
436  if(sizenew > sizeold) usenew = true;
437 
438  // samples from different events can be of different lengths - sum all
439  // that overlap.
440 
441  for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
442  if(isamp < sizenew) {
443  fC_new = (iZDC->second)[isamp];
444  }
445  else { fC_new = 0;}
446 
447  if(isamp < sizeold) {
448  fC_old = ZDC_old[isamp];
449  }
450  else { fC_old = 0;}
451 
452  // add values
453  fC_sum = fC_new + fC_old;
454 
455  if(usenew) {ZDC_bigger[isamp] = fC_sum; }
456  else { ZDC_old[isamp] = fC_sum; } // overwrite old sample, adding new info
457 
458  }
459  if(usenew) ZDC_old = ZDC_bigger; // save new, larger sized sample in "old" slot
460 
461  }
462  else {
463  if(formerID>0) {
464  // make new digi
465  ZDCdigis->push_back(ZDCDataFrame(formerID));
466 
467  // set up information to convert back
468 
469  HcalZDCDetId cell = ZDC_old.id();
470  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
471  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder); // this one is generic
472  HcalCoderDb coder (*channelCoder, *shape);
473 
474  unsigned int sizeold = ZDC_old.size();
475  for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
476  coder.fC2adc(ZDC_old,(ZDCdigis->back()), 1 ); // as per simulation, capid=1
477  }
478  }
479  //save pointers for next iteration
480  formerID = currentID;
481  ZDC_old = iZDC->second;
482  }
483 
484  iZDCchk = iZDC;
485  if((++iZDCchk) == ZDCDigiStorage_.end()) { //make sure not to lose the last one
486  // make new digi
487  ZDCdigis->push_back(ZDCDataFrame(currentID));
488 
489  // set up information to convert back
490 
491  HcalZDCDetId cell = (iZDC->second).id();
492  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
493  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder); // this one is generic
494  HcalCoderDb coder (*channelCoder, *shape);
495 
496  unsigned int sizeold = (iZDC->second).size();
497  for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
498  coder.fC2adc(ZDC_old,(ZDCdigis->back()), 1 ); // as per simulation, capid=1
499  }
500  }
501  }
502 
503 
504 
505  //done merging
506 
507  // put the collection of recunstructed hits in the event
508  LogInfo("DataMixingHcalDigiWorker") << "total # HBHE Merged digis: " << HBHEdigis->size() ;
509  LogInfo("DataMixingHcalDigiWorker") << "total # HO Merged digis: " << HOdigis->size() ;
510  LogInfo("DataMixingHcalDigiWorker") << "total # HF Merged digis: " << HFdigis->size() ;
511  LogInfo("DataMixingHcalDigiWorker") << "total # QIE10 Merged digis: " << QIE10digis->size() ;
512  LogInfo("DataMixingHcalDigiWorker") << "total # QIE11 Merged digis: " << QIE11digis->size() ;
513  LogInfo("DataMixingHcalDigiWorker") << "total # ZDC Merged digis: " << ZDCdigis->size() ;
514 
515 
516  e.put(std::move(HBHEdigis), HBHEDigiCollectionDM_);
517  e.put(std::move(HOdigis), HODigiCollectionDM_);
518  e.put(std::move(HFdigis), HFDigiCollectionDM_);
519  e.put(std::move(QIE10digis), QIE10DigiCollectionDM_);
520  e.put(std::move(QIE11digis), QIE11DigiCollectionDM_);
521  e.put(std::move(ZDCdigis), ZDCDigiCollectionDM_);
522 
523  // clear local storage after this event
524  HBHEDigiStorage_.clear();
525  HODigiStorage_.clear();
526  HFDigiStorage_.clear();
527  QIE10DigiStorage_.clear();
528  QIE11DigiStorage_.clear();
529  ZDCDigiStorage_.clear();
530 
531  }
532 
533 } //edm
#define LogDebug(id)
size
Write out results.
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
edm::EDGetTokenT< HBHEDigiCollection > HBHEDigiPToken_
edm::EDGetTokenT< HODigiCollection > HODigiToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
edm::EDGetTokenT< ZDCDigiCollection > ZDCDigiToken_
void setSize(unsigned int size)
Definition: CaloSamples.h:55
edm::EDGetTokenT< HODigiCollection > HODigiPToken_
EventID const & id() const
std::vector< ZDCDataFrame >::const_iterator const_iterator
double pedestal(int fCapId) const
get pedestal for capid=0..3
void addHcalPileups(const int bcr, const edm::EventPrincipal *, unsigned int EventId, const edm::EventSetup &ES, ModuleCallingContext const *)
virtual void fC2adc(const CaloSamples &clf, HBHEDataFrame &df, int fCapIdOffset) const
Definition: HcalCoderDb.cc:76
edm::EDGetTokenT< ZDCDigiCollection > ZDCDigiPToken_
edm::EDGetTokenT< HFDigiCollection > HFDigiPToken_
void putHcal(edm::Event &e, const edm::EventSetup &ES)
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const
Definition: HcalCoderDb.cc:68
const_iterator end() const
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
int size() const
get the size
Definition: CaloSamples.h:24
const T & get() const
Definition: EventSetup.h:56
const HcalQIECoder * getHcalCoder(const HcalGenericDetId &fId) const
edm::EDGetTokenT< HFDigiCollection > HFDigiToken_
const HcalQIEShape * getHcalShape(const HcalGenericDetId &fId) const
std::string const & label() const
Definition: InputTag.h:36
edm::EDGetTokenT< QIE11DigiCollection > QIE11DigiToken_
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
size_type size() const
DetId id() const
get the (generic) id
Definition: CaloSamples.h:21
edm::EDGetTokenT< QIE11DigiCollection > QIE11DigiPToken_
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
void addHcalSignals(const edm::Event &e, const edm::EventSetup &ES)
edm::EDGetTokenT< QIE10DigiCollection > QIE10DigiToken_
edm::EDGetTokenT< QIE10DigiCollection > QIE10DigiPToken_
def move(src, dest)
Definition: eostools.py:510
const_iterator begin() const
edm::EDGetTokenT< HBHEDigiCollection > HBHEDigiToken_