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