CMS 3D CMS Logo

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