CMS 3D CMS Logo

DataMixingEMDigiWorker.cc
Go to the documentation of this file.
1 // File: DataMixingEMDigiWorker.cc
2 // Description: see DataMixingEMDigiWorker.h
3 // Author: Mike Hildreth, University of Notre Dame
4 //
5 //--------------------------------------------
6 
16 #include <cmath>
17 #include <map>
18 #include <memory>
19 
20 //
21 //
22 #include "DataMixingEMDigiWorker.h"
23 
24 using namespace std;
25 
26 namespace edm {
27 
28  // Virtual constructor
29 
30  DataMixingEMDigiWorker::DataMixingEMDigiWorker() {}
31 
32  // Constructor
33  DataMixingEMDigiWorker::DataMixingEMDigiWorker(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
34  : label_(ps.getParameter<std::string>("Label"))
35 
36  {
37  // get the subdetector names
38  // this->getSubdetectorNames(); //something like this may be useful to
39  // check what we are supposed to do...
40 
41  // declare the products to produce, retrieve
42 
43  EBProducerSig_ = ps.getParameter<edm::InputTag>("EBdigiProducerSig");
44  EEProducerSig_ = ps.getParameter<edm::InputTag>("EEdigiProducerSig");
45  ESProducerSig_ = ps.getParameter<edm::InputTag>("ESdigiProducerSig");
46 
50 
51  EBPileInputTag_ = ps.getParameter<edm::InputTag>("EBPileInputTag");
52  EEPileInputTag_ = ps.getParameter<edm::InputTag>("EEPileInputTag");
53  ESPileInputTag_ = ps.getParameter<edm::InputTag>("ESPileInputTag");
54 
58 
59  EBDigiCollectionDM_ = ps.getParameter<std::string>("EBDigiCollectionDM");
60  EEDigiCollectionDM_ = ps.getParameter<std::string>("EEDigiCollectionDM");
61  ESDigiCollectionDM_ = ps.getParameter<std::string>("ESDigiCollectionDM");
62  }
63 
64  // Virtual destructor needed.
66 
68  // fill in maps of hits
69 
70  LogInfo("DataMixingEMDigiWorker") << "===============> adding MC signals for " << e.id();
71 
72  // EB first
73 
74  Handle<EBDigiCollection> pEBDigis;
75 
76  const EBDigiCollection *EBDigis = nullptr;
77 
78  if (e.getByToken(EBDigiToken_, pEBDigis)) {
79  EBDigis = pEBDigis.product(); // get a ptr to the product
80  LogDebug("DataMixingEMDigiWorker") << "total # EB digis: " << EBDigis->size();
81  }
82 
83  if (EBDigis) {
84  // loop over digis, storing them in a map so we can add pileup later
85 
86  for (EBDigiCollection::const_iterator it = EBDigis->begin(); it != EBDigis->end(); ++it) {
87  EBDigiStorage_.insert(EBDigiMap::value_type((it->id()), *it));
88 #ifdef DEBUG
89  // Commented out because this does not compile anymore
90  // LogDebug("DataMixingEMDigiWorker") << "processed EBDigi with rawId: "
91  // << it->id().rawId() << "\n"
92  // << " digi energy: " << it->energy();
93 #endif
94  }
95  }
96 
97  // EE next
98 
99  Handle<EEDigiCollection> pEEDigis;
100 
101  const EEDigiCollection *EEDigis = nullptr;
102 
103  if (e.getByToken(EEDigiToken_, pEEDigis)) {
104  EEDigis = pEEDigis.product(); // get a ptr to the product
105  LogDebug("DataMixingEMDigiWorker") << "total # EE digis: " << EEDigis->size();
106  }
107 
108  if (EEDigis) {
109  // loop over digis, storing them in a map so we can add pileup later
110  for (EEDigiCollection::const_iterator it = EEDigis->begin(); it != EEDigis->end(); ++it) {
111  EEDigiStorage_.insert(EEDigiMap::value_type((it->id()), *it));
112 #ifdef DEBUG
113  // Commented out because this does not compile anymore
114  // LogDebug("DataMixingEMDigiWorker") << "processed EEDigi with rawId: "
115  // << it->id().rawId() << "\n"
116  // << " digi energy: " << it->energy();
117 #endif
118  }
119  }
120  // ES next
121 
122  Handle<ESDigiCollection> pESDigis;
123 
124  const ESDigiCollection *ESDigis = nullptr;
125 
126  if (e.getByToken(ESDigiToken_, pESDigis)) {
127  ESDigis = pESDigis.product(); // get a ptr to the product
128 #ifdef DEBUG
129  LogDebug("DataMixingEMDigiWorker") << "total # ES digis: " << ESDigis->size();
130 #endif
131  }
132 
133  if (ESDigis) {
134  // loop over digis, storing them in a map so we can add pileup later
135  for (ESDigiCollection::const_iterator it = ESDigis->begin(); it != ESDigis->end(); ++it) {
136  ESDigiStorage_.insert(ESDigiMap::value_type((it->id()), *it));
137 
138 #ifdef DEBUG
139  // Commented out because this does not compile anymore
140  // LogDebug("DataMixingEMDigiWorker") << "processed ESDigi with rawId: "
141  // << it->id().rawId() << "\n"
142  // << " digi energy: " << it->energy();
143 #endif
144  }
145  }
146 
147  } // end of addEMSignals
148 
150  const EventPrincipal *ep,
151  unsigned int eventNr,
152  const edm::EventSetup &ES,
153  ModuleCallingContext const *mcc) {
154  LogInfo("DataMixingEMDigiWorker") << "\n===============> adding pileups from event " << ep->id()
155  << " for bunchcrossing " << bcr;
156 
157  // fill in maps of hits; same code as addSignals, except now applied to the
158  // pileup events
159 
160  // EB first
161 
162  std::shared_ptr<Wrapper<EBDigiCollection> const> EBDigisPTR =
163  getProductByTag<EBDigiCollection>(*ep, EBPileInputTag_, mcc);
164 
165  if (EBDigisPTR) {
166  const EBDigiCollection *EBDigis = const_cast<EBDigiCollection *>(EBDigisPTR->product());
167 
168  LogDebug("DataMixingEMDigiWorker") << "total # EB digis: " << EBDigis->size();
169 
170  // loop over digis, adding these to the existing maps
171  for (EBDigiCollection::const_iterator it = EBDigis->begin(); it != EBDigis->end(); ++it) {
172  EBDigiStorage_.insert(EBDigiMap::value_type((it->id()), *it));
173 
174 #ifdef DEBUG
175  // Commented out because this does not compile anymore
176  // LogDebug("DataMixingEMDigiWorker") << "processed EBDigi with rawId: "
177  // << it->id().rawId() << "\n"
178  // << " digi energy: " << it->energy();
179 #endif
180  }
181  }
182 
183  // EE Next
184 
185  std::shared_ptr<Wrapper<EEDigiCollection> const> EEDigisPTR =
186  getProductByTag<EEDigiCollection>(*ep, EEPileInputTag_, mcc);
187 
188  if (EEDigisPTR) {
189  const EEDigiCollection *EEDigis = const_cast<EEDigiCollection *>(EEDigisPTR->product());
190 
191  LogDebug("DataMixingEMDigiWorker") << "total # EE digis: " << EEDigis->size();
192 
193  for (EEDigiCollection::const_iterator it = EEDigis->begin(); it != EEDigis->end(); ++it) {
194  EEDigiStorage_.insert(EEDigiMap::value_type((it->id()), *it));
195 
196 #ifdef DEBUG
197  // Commented out because this does not compile anymore
198  // LogDebug("DataMixingEMDigiWorker") << "processed EEDigi with rawId: "
199  // << it->id().rawId() << "\n"
200  // << " digi energy: " << it->energy();
201 #endif
202  }
203  }
204  // ES Next
205 
206  std::shared_ptr<Wrapper<ESDigiCollection> const> ESDigisPTR =
207  getProductByTag<ESDigiCollection>(*ep, ESPileInputTag_, mcc);
208 
209  if (ESDigisPTR) {
210  const ESDigiCollection *ESDigis = const_cast<ESDigiCollection *>(ESDigisPTR->product());
211 
212  LogDebug("DataMixingEMDigiWorker") << "total # ES digis: " << ESDigis->size();
213 
214  for (ESDigiCollection::const_iterator it = ESDigis->begin(); it != ESDigis->end(); ++it) {
215  ESDigiStorage_.insert(ESDigiMap::value_type((it->id()), *it));
216 
217 #ifdef DEBUG
218  // Commented out because this does not compile anymore
219  // LogDebug("DataMixingEMDigiWorker") << "processed ESDigi with rawId: "
220  // << it->id().rawId() << "\n"
221  // << " digi energy: " << it->energy();
222 #endif
223  }
224  }
225  }
226 
228  // collection of digis to put in the event
229  std::unique_ptr<EBDigiCollection> EBdigis(new EBDigiCollection);
230  std::unique_ptr<EEDigiCollection> EEdigis(new EEDigiCollection);
231  std::unique_ptr<ESDigiCollection> ESdigis(new ESDigiCollection);
232 
233  // loop over the maps we have, re-making individual hits or digis if
234  // necessary.
235  DetId formerID = 0;
236  DetId currentID;
237 
238  EBDataFrame EB_old;
239 
240  int gain_new = 0;
241  int gain_old = 0;
242  int gain_consensus = 0;
243  int adc_new;
244  int adc_old;
245  int adc_sum;
246  uint16_t data;
247 
248  // EB first...
249 
250  EBDigiMap::const_iterator iEBchk;
251 
252  for (EBDigiMap::const_iterator iEB = EBDigiStorage_.begin(); iEB != EBDigiStorage_.end(); iEB++) {
253  currentID = iEB->first;
254 
255  if (currentID == formerID) { // we have to add these digis together
256  /*
257  cout<< " Adding signals " << EBDetId(currentID).ieta() << " "
258  << EBDetId(currentID).iphi() << std::endl;
259 
260  cout << 1 << " " ;
261  for (int i=0; i<10;++i) std::cout << EB_old[i].adc()<<
262  "["<<EB_old[i].gainId()<< "] " ; std::cout << std::endl;
263 
264  cout << 2 << " " ;
265  for (int i=0; i<10;++i) std::cout << (iEB->second)[i].adc()<<
266  "["<<(iEB->second)[i].gainId()<< "] " ; std::cout << std::endl;
267  */
268  // loop over digi samples in each DataFrame
269  unsigned int sizenew = (iEB->second).size();
270  unsigned int sizeold = EB_old.size();
271 
272  unsigned int max_samp = std::max(sizenew, sizeold);
273 
274  // samples from different events can be of different lengths - sum all
275  // that overlap.
276  // check to see if gains match - if not, scale smaller cell down.
277 
278  int sw_gain_consensus = 0;
279 
280  for (unsigned int isamp = 0; isamp < max_samp; isamp++) {
281  if (isamp < sizenew) {
282  gain_new = (iEB->second)[isamp].gainId();
283  adc_new = (iEB->second)[isamp].adc();
284  } else {
285  adc_new = 0;
286  }
287 
288  if (isamp < sizeold) {
289  gain_old = EB_old[isamp].gainId();
290  adc_old = EB_old[isamp].adc();
291  } else {
292  adc_old = 0;
293  }
294 
295  const std::vector<float> pedeStals = GetPedestals(ES, currentID);
296  const std::vector<float> gainRatios = GetGainRatios(ES, currentID);
297 
298  if (adc_new > 0 && adc_old > 0) {
299  if (gain_old == gain_new) { // we're happy - easy case
300  gain_consensus = gain_old;
301  } else { // lower gain sample has more energy
302 
303  if (gain_old < gain_new) { // old has higher gain than new, scale to lower gain
304 
305  float ratio = gainRatios[gain_new - 1] / gainRatios[gain_old - 1];
306  adc_old = (int)round((adc_old - pedeStals[gain_old - 1]) / ratio + pedeStals[gain_new - 1]);
307  gain_consensus = gain_new;
308  } else { // scale to old (lower) gain
309  float ratio = gainRatios[gain_old - 1] / gainRatios[gain_new - 1];
310  adc_new = (int)round((adc_new - pedeStals[gain_new - 1]) / ratio + pedeStals[gain_old - 1]);
311  gain_consensus = gain_old;
312  }
313  }
314  }
315 
316  // add values, but don't count pedestals twice
317  adc_sum = adc_new + adc_old - (int)round(pedeStals[gain_consensus - 1]);
318 
319  // if we are now saturating that gain, switch to the next
320  if (adc_sum > 4096) {
321  if (gain_consensus < 3) {
322  double ratio = gainRatios[gain_consensus] / gainRatios[gain_consensus - 1];
323  adc_sum = (int)round((adc_sum - pedeStals[gain_consensus - 1]) / ratio + pedeStals[gain_consensus]);
324  sw_gain_consensus = ++gain_consensus;
325  } else
326  adc_sum = 4096;
327  }
328 
329  // furthermore, make sure we don't decrease our gain once we've switched
330  // up in case go back
331  if (gain_consensus < sw_gain_consensus) {
332  double ratio = gainRatios[sw_gain_consensus - 1] / gainRatios[gain_consensus - 1];
333  adc_sum = (int)round((adc_sum - pedeStals[gain_consensus - 1]) / ratio + pedeStals[sw_gain_consensus - 1]);
334  gain_consensus = sw_gain_consensus;
335  }
336 
337  EcalMGPASample sample(adc_sum, gain_consensus);
338  EB_old.setSample(isamp,
339  sample); // overwrite old sample, adding new info
340  } // for sample
341 
342  } // if current = former
343  else {
344  if (formerID > 0) {
345  EBdigis->push_back(formerID, EB_old.frame().begin());
346  }
347  // save pointers for next iteration
348  formerID = currentID;
349  EB_old = iEB->second;
350  }
351 
352  iEBchk = iEB;
353  if ((++iEBchk) == EBDigiStorage_.end()) { // make sure not to lose the last one
354  EBdigis->push_back(currentID, (iEB->second).frame().begin());
355  }
356  }
357 
358  // EE next...
359 
360  formerID = 0;
361  EEDataFrame EE_old;
362 
363  EEDigiMap::const_iterator iEEchk;
364 
365  for (EEDigiMap::const_iterator iEE = EEDigiStorage_.begin(); iEE != EEDigiStorage_.end(); iEE++) {
366  currentID = iEE->first;
367 
368  if (currentID == formerID) { // we have to add these digis together
369 
370  // loop over digi samples in each DataFrame
371  unsigned int sizenew = (iEE->second).size();
372  unsigned int sizeold = EE_old.size();
373 
374  unsigned int max_samp = std::max(sizenew, sizeold);
375 
376  // samples from different events can be of different lengths - sum all
377  // that overlap.
378  // check to see if gains match - if not, scale smaller cell down.
379 
380  for (unsigned int isamp = 0; isamp < max_samp; isamp++) {
381  if (isamp < sizenew) {
382  gain_new = (iEE->second)[isamp].gainId();
383  adc_new = (iEE->second)[isamp].adc();
384  } else {
385  adc_new = 0;
386  }
387 
388  if (isamp < sizeold) {
389  gain_old = EE_old[isamp].gainId();
390  adc_old = EE_old[isamp].adc();
391  } else {
392  adc_old = 0;
393  }
394 
395  const std::vector<float> pedeStals = GetPedestals(ES, currentID);
396  const std::vector<float> gainRatios = GetGainRatios(ES, currentID);
397 
398  if (adc_new > 0 && adc_old > 0) {
399  if (gain_old == gain_new) { // we're happy - easy case
400  gain_consensus = gain_old;
401  } else { // lower gain sample has more energy
402 
403  if (gain_old < gain_new) { // old has higher gain than new, scale to lower gain
404 
405  float ratio = gainRatios[gain_new - 1] / gainRatios[gain_old - 1];
406  adc_old = (int)round((adc_old - pedeStals[gain_old - 1]) / ratio + pedeStals[gain_new - 1]);
407  gain_consensus = gain_new;
408  } else { // scale to old (lower) gain
409  float ratio = gainRatios[gain_old - 1] / gainRatios[gain_new - 1];
410  adc_new = (int)round((adc_new - pedeStals[gain_new - 1]) / ratio + pedeStals[gain_old - 1]);
411  gain_consensus = gain_old;
412  }
413  }
414  }
415 
416  // add values, but don't count pedestals twice
417  adc_sum = adc_new + adc_old - (int)round(pedeStals[gain_consensus - 1]);
418 
419  // if the sum saturates this gain, switch
420  if (adc_sum > 4096) {
421  if (gain_consensus < 3) {
422  double ratio = gainRatios[gain_consensus] / gainRatios[gain_consensus - 1];
423  adc_sum = (int)round((adc_sum - pedeStals[gain_consensus - 1]) / ratio + pedeStals[gain_consensus]);
424  ++gain_consensus;
425  } else
426  adc_sum = 4096;
427  }
428 
429  EcalMGPASample sample(adc_sum, gain_consensus);
430  EE_old.setSample(isamp, sample);
431  }
432 
433  } else {
434  if (formerID > 0) {
435  EEdigis->push_back(formerID, EE_old.frame().begin());
436  }
437  // save pointers for next iteration
438  formerID = currentID;
439  EE_old = iEE->second;
440  }
441 
442  iEEchk = iEE;
443  if ((++iEEchk) == EEDigiStorage_.end()) { // make sure not to lose the last one
444  EEdigis->push_back(currentID, (iEE->second).frame().begin());
445  }
446  }
447 
448  // ES next...
449 
450  formerID = 0;
451  ESDataFrame ES_old;
452 
453  ESDigiMap::const_iterator iESchk;
454 
455  for (ESDigiMap::const_iterator iES = ESDigiStorage_.begin(); iES != ESDigiStorage_.end(); iES++) {
456  currentID = iES->first;
457 
458  if (currentID == formerID) { // we have to add these digis together
459 
460  // loop over digi samples in each DataFrame
461  unsigned int sizenew = (iES->second).size();
462  unsigned int sizeold = ES_old.size();
463  uint16_t rawdat = 0;
464  unsigned int max_samp = std::max(sizenew, sizeold);
465 
466  // samples from different events can be of different lengths - sum all
467  // that overlap.
468  // check to see if gains match - if not, scale smaller cell down.
469 
470  for (unsigned int isamp = 0; isamp < max_samp; isamp++) {
471  if (isamp < sizenew) {
472  adc_new = (iES->second)[isamp].adc();
473  rawdat = (iES->second)[isamp].raw();
474  } else {
475  adc_new = 0;
476  }
477 
478  if (isamp < sizeold) {
479  adc_old = ES_old[isamp].adc();
480  rawdat = ES_old[isamp].raw();
481  } else {
482  adc_old = 0;
483  }
484 
485  // add values
486  adc_sum = adc_new + adc_old;
487  // make data word of gain, rawdata
488  adc_sum = std::min(adc_sum, 4095); // first 12 bits of (uint)
489  data = adc_sum + (rawdat & 0xF000); // data is 14 bit word with gain as MSBs
490  ES_old.setSample(isamp, data);
491  }
492 
493  } else {
494  if (formerID > 0) {
495  ESdigis->push_back(ES_old);
496  }
497  // save pointers for next iteration
498  formerID = currentID;
499  ES_old = iES->second;
500  }
501 
502  iESchk = iES;
503  if ((++iESchk) == ESDigiStorage_.end()) { // make sure not to lose the last one
504  ESdigis->push_back(iES->second);
505  // ESDataFrame df( (*ESdigis)->back() );
506  // for(int isamp=0; isamp<(iES->second).size(); isamp++) {
507  // df.setSample(isamp,(iES->second).data[isamp]);
508  // }
509  }
510  }
511 
512  // done merging
513 
514  // put the collection of reconstructed hits in the event
515  LogInfo("DataMixingEMDigiWorker") << "total # EB Merged digis: " << EBdigis->size();
516  LogInfo("DataMixingEMDigiWorker") << "total # EE Merged digis: " << EEdigis->size();
517  LogInfo("DataMixingEMDigiWorker") << "total # ES Merged digis: " << ESdigis->size();
518 
519  e.put(std::move(EBdigis), EBDigiCollectionDM_);
520  e.put(std::move(EEdigis), EEDigiCollectionDM_);
521  e.put(std::move(ESdigis), ESDigiCollectionDM_);
522 
523  // clear local storage after this event
524 
525  EBDigiStorage_.clear();
526  EEDigiStorage_.clear();
527  ESDigiStorage_.clear();
528  }
529  const std::vector<float> DataMixingEMDigiWorker::GetPedestals(const edm::EventSetup &ES, const DetId &detid) {
530  std::vector<float> pedeStals(3);
531 
532  // get pedestals
534  ES.get<EcalPedestalsRcd>().get(pedHandle);
535 
536  const EcalPedestalsMap &pedMap = pedHandle.product()->getMap(); // map of pedestals
537  EcalPedestalsMapIterator pedIter; // pedestal iterator
538  EcalPedestals::Item aped; // pedestal object for a single xtal
539 
540  pedIter = pedMap.find(detid);
541  if (pedIter != pedMap.end()) {
542  aped = (*pedIter);
543  pedeStals[0] = aped.mean_x12;
544  pedeStals[1] = aped.mean_x6;
545  pedeStals[2] = aped.mean_x1;
546  } else {
547  edm::LogError("DataMixingMissingInput") << "Cannot find pedestals";
548  pedeStals[0] = 0;
549  pedeStals[1] = 0;
550  pedeStals[2] = 0;
551  }
552 
553  return pedeStals;
554  }
555 
556  const std::vector<float> DataMixingEMDigiWorker::GetGainRatios(const edm::EventSetup &ES, const DetId &detid) {
557  std::vector<float> gainRatios(3);
558  // get gain ratios
560  ES.get<EcalGainRatiosRcd>().get(grHandle);
561  EcalMGPAGainRatio theRatio = (*grHandle)[detid];
562 
563  gainRatios[0] = 1.;
564  gainRatios[1] = theRatio.gain12Over6();
565  gainRatios[2] = theRatio.gain6Over1() * theRatio.gain12Over6();
566 
567  return gainRatios;
568  }
569 
570 } // namespace 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:125
const std::vector< float > GetGainRatios(const edm::EventSetup &ES, const DetId &detid)
retrieve gain ratios for that detid [0]=g12, [1]=g6, [2]=g12
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
edm::EDGetTokenT< EEDigiCollection > EEDigiPileToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< ESDigiCollection > ESDigiPileToken_
const self & getMap() const
EventID const & id() const
const std::vector< float > GetPedestals(const edm::EventSetup &ES, const DetId &detid)
retrieve pedestals for that detid [0]=g12, [1]=g6, [2]=g12
const_iterator begin() const
int size() const
Definition: ESDataFrame.h:23
int size() const
Definition: EcalDataFrame.h:26
void addEMPileups(const int bcr, const edm::EventPrincipal *, unsigned int EventId, const edm::EventSetup &ES, ModuleCallingContext const *)
void putEM(edm::Event &e, const edm::EventSetup &ES)
EcalPedestalsMap::const_iterator EcalPedestalsMapIterator
Definition: EcalPedestals.h:52
float gain6Over1() const
T min(T a, T b)
Definition: MathUtil.h:58
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
void setSample(int i, const ESSample &sam)
Definition: ESDataFrame.h:30
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:74
constexpr iterator begin()
Definition: DataFrame.h:48
edm::DataFrame const & frame() const
Definition: EcalDataFrame.h:50
float gain12Over6() const
const_iterator end() const
edm::EventID id() const
Definition: EventBase.h:59
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void setSample(int i, EcalMGPASample sam)
Definition: EcalDataFrame.h:43
T get() const
Definition: EventSetup.h:71
const_iterator find(uint32_t rawId) const
const_iterator end() const
void addEMSignals(const edm::Event &e, const edm::EventSetup &ES)
T const * product() const
Definition: ESHandle.h:86
edm::EDGetTokenT< EBDigiCollection > EBDigiToken_
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< ESDigiCollection > ESDigiToken_
edm::EDGetTokenT< EEDigiCollection > EEDigiToken_
edm::EDGetTokenT< EBDigiCollection > EBDigiPileToken_