CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/SimGeneral/DataMixingModule/plugins/DataMixingEMDigiWorker.cc

Go to the documentation of this file.
00001 // File: DataMixingEMDigiWorker.cc
00002 // Description:  see DataMixingEMDigiWorker.h
00003 // Author:  Mike Hildreth, University of Notre Dame
00004 //
00005 //--------------------------------------------
00006 
00007 #include <map>
00008 #include <cmath>
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "FWCore/Utilities/interface/EDMException.h"
00011 #include "FWCore/Framework/interface/ConstProductRegistry.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "DataFormats/Common/interface/Handle.h"
00014 #include "DataFormats/Provenance/interface/Provenance.h"
00015 #include "DataFormats/Provenance/interface/BranchDescription.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "CondFormats/EcalObjects/interface/EcalGainRatios.h"
00019 #include "CondFormats/DataRecord/interface/EcalGainRatiosRcd.h"
00020 #include "CondFormats/EcalObjects/interface/EcalPedestals.h"
00021 #include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
00022 #include "CondFormats/EcalObjects/interface/EcalMGPAGainRatio.h"
00023 
00024 //
00025 //
00026 #include "DataMixingEMDigiWorker.h"
00027 
00028 
00029 using namespace std;
00030 
00031 namespace edm
00032 {
00033 
00034   // Virtual constructor
00035 
00036   DataMixingEMDigiWorker::DataMixingEMDigiWorker() { }
00037 
00038   // Constructor 
00039   DataMixingEMDigiWorker::DataMixingEMDigiWorker(const edm::ParameterSet& ps) : 
00040                                                             label_(ps.getParameter<std::string>("Label"))
00041 
00042   {                                                         
00043 
00044     // get the subdetector names
00045     //    this->getSubdetectorNames();  //something like this may be useful to check what we are supposed to do...
00046 
00047     // declare the products to produce, retrieve
00048 
00049     EBProducerSig_ = ps.getParameter<edm::InputTag>("EBdigiProducerSig");
00050     EEProducerSig_ = ps.getParameter<edm::InputTag>("EEdigiProducerSig");
00051     ESProducerSig_ = ps.getParameter<edm::InputTag>("ESdigiProducerSig");
00052 
00053     EBPileInputTag_ = ps.getParameter<edm::InputTag>("EBPileInputTag");
00054     EEPileInputTag_ = ps.getParameter<edm::InputTag>("EEPileInputTag");
00055     ESPileInputTag_ = ps.getParameter<edm::InputTag>("ESPileInputTag");
00056 
00057     EBDigiCollectionDM_        = ps.getParameter<std::string>("EBDigiCollectionDM");
00058     EEDigiCollectionDM_        = ps.getParameter<std::string>("EEDigiCollectionDM");
00059     ESDigiCollectionDM_        = ps.getParameter<std::string>("ESDigiCollectionDM");
00060    //   nMaxPrintout_            = ps.getUntrackedParameter<int>("nMaxPrintout",10);
00061 
00062    //EBalgo_ = new EcalDigiSimpleAlgo();
00063    //EEalgo_ = new EcalDigiSimpleAlgo();
00064 
00065    // don't think I can "produce" in a sub-class...
00066 
00067    //produces< EBDigiCollection >(EBDigiCollectionDM_);
00068    //produces< EEDigiCollection >(EEDigiCollectionDM_);
00069 
00070   }
00071                
00072 
00073   // Virtual destructor needed.
00074   DataMixingEMDigiWorker::~DataMixingEMDigiWorker() { 
00075   }  
00076 
00077   void DataMixingEMDigiWorker::addEMSignals(const edm::Event &e,const edm::EventSetup& ES) { 
00078     // fill in maps of hits
00079 
00080     LogInfo("DataMixingEMDigiWorker")<<"===============> adding MC signals for "<<e.id();
00081 
00082     // EB first
00083 
00084    Handle< EBDigiCollection > pEBDigis;
00085 
00086    const EBDigiCollection*  EBDigis = 0;
00087 
00088    if(e.getByLabel(EBProducerSig_, pEBDigis) ){
00089      EBDigis = pEBDigis.product(); // get a ptr to the product
00090      LogDebug("DataMixingEMDigiWorker") << "total # EB digis: " << EBDigis->size();
00091    }
00092    else { std::cout << "NO EBDigis! " << EBProducerSig_.label() << " " << EBdigiCollectionSig_.label() << std::endl;}
00093  
00094    if (EBDigis)
00095      {
00096        // loop over digis, storing them in a map so we can add pileup later
00097 
00098 
00099        for(EBDigiCollection::const_iterator it  = EBDigis->begin();     
00100            it != EBDigis->end(); ++it) {
00101 
00102          EBDigiStorage_.insert(EBDigiMap::value_type( ( it->id() ), *it ));
00103 #ifdef DEBUG     
00104          LogDebug("DataMixingEMDigiWorker") << "processed EBDigi with rawId: "
00105                                       << it->id().rawId() << "\n"
00106                                       << " digi energy: " << it->energy();
00107 #endif
00108        }
00109      }
00110 
00111    // EE next
00112 
00113    Handle< EEDigiCollection > pEEDigis;
00114 
00115    const EEDigiCollection*  EEDigis = 0;
00116 
00117    
00118    if(e.getByLabel(EEProducerSig_, pEEDigis) ){
00119      EEDigis = pEEDigis.product(); // get a ptr to the product
00120      LogDebug("DataMixingEMDigiWorker") << "total # EE digis: " << EEDigis->size();
00121    } 
00122    
00123  
00124    if (EEDigis)
00125      {
00126        // loop over digis, storing them in a map so we can add pileup later
00127        for(EEDigiCollection::const_iterator it  = EEDigis->begin();     
00128            it != EEDigis->end(); ++it) {
00129 
00130          EEDigiStorage_.insert(EEDigiMap::value_type( ( it->id() ), *it ));
00131 #ifdef DEBUG     
00132          LogDebug("DataMixingEMDigiWorker") << "processed EEDigi with rawId: "
00133                                       << it->id().rawId() << "\n"
00134                                       << " digi energy: " << it->energy();
00135 #endif
00136 
00137        }
00138      }
00139    // ES next
00140 
00141    Handle< ESDigiCollection > pESDigis;
00142 
00143    const ESDigiCollection*  ESDigis = 0;
00144 
00145    
00146    if(e.getByLabel( ESProducerSig_, pESDigis) ){
00147      ESDigis = pESDigis.product(); // get a ptr to the product
00148 #ifdef DEBUG
00149      LogDebug("DataMixingEMDigiWorker") << "total # ES digis: " << ESDigis->size();
00150 #endif
00151    } 
00152 
00153  
00154    if (ESDigis)
00155      {
00156 
00157        // loop over digis, storing them in a map so we can add pileup later
00158        for(ESDigiCollection::const_iterator it  = ESDigis->begin();     
00159            it != ESDigis->end(); ++it) {
00160 
00161          ESDigiStorage_.insert(ESDigiMap::value_type( ( it->id() ), *it ));
00162          
00163 #ifdef DEBUG     
00164          LogDebug("DataMixingEMDigiWorker") << "processed ESDigi with rawId: "
00165                                       << it->id().rawId() << "\n"
00166                                       << " digi energy: " << it->energy();
00167 #endif
00168 
00169        }
00170      }
00171     
00172   } // end of addEMSignals
00173 
00174   void DataMixingEMDigiWorker::addEMPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr, const edm::EventSetup& ES) {
00175   
00176     LogInfo("DataMixingEMDigiWorker") <<"\n===============> adding pileups from event  "<<ep->id()<<" for bunchcrossing "<<bcr;
00177 
00178     // fill in maps of hits; same code as addSignals, except now applied to the pileup events
00179 
00180     // EB first
00181 
00182     boost::shared_ptr<Wrapper<EBDigiCollection>  const> EBDigisPTR = 
00183           getProductByTag<EBDigiCollection>(*ep, EBPileInputTag_ );
00184  
00185    if(EBDigisPTR ) {
00186 
00187      const EBDigiCollection*  EBDigis = const_cast< EBDigiCollection * >(EBDigisPTR->product());
00188 
00189      LogDebug("DataMixingEMDigiWorker") << "total # EB digis: " << EBDigis->size();
00190 
00191        // loop over digis, adding these to the existing maps
00192      for(EBDigiCollection::const_iterator it  = EBDigis->begin();
00193          it != EBDigis->end(); ++it) {
00194 
00195          EBDigiStorage_.insert(EBDigiMap::value_type( (it->id()), *it ));
00196          
00197 #ifdef DEBUG     
00198          LogDebug("DataMixingEMDigiWorker") << "processed EBDigi with rawId: "
00199                       << it->id().rawId() << "\n"
00200                       << " digi energy: " << it->energy();
00201 #endif
00202       }
00203     }
00204 
00205     // EE Next
00206 
00207     boost::shared_ptr<Wrapper<EEDigiCollection>  const> EEDigisPTR =
00208           getProductByTag<EEDigiCollection>(*ep, EEPileInputTag_ );
00209 
00210     if(EEDigisPTR ) {
00211 
00212      const EEDigiCollection*  EEDigis = const_cast< EEDigiCollection * >(EEDigisPTR->product()); 
00213 
00214      LogDebug("DataMixingEMDigiWorker") << "total # EE digis: " << EEDigis->size();
00215 
00216      for(EEDigiCollection::const_iterator it  = EEDigis->begin();
00217          it != EEDigis->end(); ++it) {
00218 
00219        EEDigiStorage_.insert(EEDigiMap::value_type( (it->id()), *it ));
00220          
00221 #ifdef DEBUG     
00222        LogDebug("DataMixingEMDigiWorker") << "processed EEDigi with rawId: "
00223                                       << it->id().rawId() << "\n"
00224                                       << " digi energy: " << it->energy();
00225 #endif
00226      }
00227    }
00228     // ES Next
00229 
00230     boost::shared_ptr<Wrapper<ESDigiCollection>  const> ESDigisPTR =
00231       getProductByTag<ESDigiCollection>(*ep, ESPileInputTag_ );
00232 
00233     if(ESDigisPTR ) {
00234 
00235       const ESDigiCollection*  ESDigis = const_cast< ESDigiCollection * >(ESDigisPTR->product());
00236 
00237       LogDebug("DataMixingEMDigiWorker") << "total # ES digis: " << ESDigis->size();
00238 
00239       for(ESDigiCollection::const_iterator it  = ESDigis->begin();
00240           it != ESDigis->end(); ++it) {
00241 
00242         ESDigiStorage_.insert(ESDigiMap::value_type( (it->id()), *it ));
00243 
00244 #ifdef DEBUG
00245         LogDebug("DataMixingEMDigiWorker") << "processed ESDigi with rawId: "
00246                                            << it->id().rawId() << "\n"
00247                                            << " digi energy: " << it->energy();
00248 #endif
00249       }
00250     }
00251  
00252   
00253  
00254   }
00255  
00256   void DataMixingEMDigiWorker::putEM(edm::Event &e, const edm::EventSetup& ES) {
00257 
00258     // collection of digis to put in the event
00259     std::auto_ptr< EBDigiCollection > EBdigis( new EBDigiCollection );
00260     std::auto_ptr< EEDigiCollection > EEdigis( new EEDigiCollection );
00261     std::auto_ptr< ESDigiCollection > ESdigis( new ESDigiCollection );
00262 
00263 
00264     // loop over the maps we have, re-making individual hits or digis if necessary.
00265     DetId formerID = 0;
00266     DetId currentID;
00267 
00268     EBDataFrame EB_old;
00269 
00270     int gain_new = 0;
00271     int gain_old = 0;
00272     int gain_consensus = 0;
00273     int adc_new;
00274     int adc_old;
00275     int adc_sum;
00276     uint16_t data;
00277 
00278     // EB first...
00279 
00280     EBDigiMap::const_iterator iEBchk;
00281 
00282     
00283 
00284     for(EBDigiMap::const_iterator iEB  = EBDigiStorage_.begin();
00285         iEB != EBDigiStorage_.end(); iEB++) {
00286 
00287       currentID = iEB->first; 
00288 
00289       if (currentID == formerID) { // we have to add these digis together
00290         /*      
00291         cout<< " Adding signals " << EBDetId(currentID).ieta() << " " 
00292                                   << EBDetId(currentID).iphi() << std::endl;
00293 
00294         cout << 1 << " " ; 
00295         for (int i=0; i<10;++i)  std::cout << EB_old[i].adc()<< "["<<EB_old[i].gainId()<< "] " ; std::cout << std::endl;
00296  
00297         cout << 2 << " " ; 
00298         for (int i=0; i<10;++i)  std::cout << (iEB->second)[i].adc()<< "["<<(iEB->second)[i].gainId()<< "] " ; std::cout << std::endl;
00299         */
00300         //loop over digi samples in each DataFrame
00301         unsigned int sizenew = (iEB->second).size();
00302         unsigned int sizeold = EB_old.size();
00303 
00304         unsigned int max_samp = std::max(sizenew, sizeold);
00305 
00306         // samples from different events can be of different lengths - sum all
00307         // that overlap.
00308         // check to see if gains match - if not, scale smaller cell down.
00309 
00310         int sw_gain_consensus=0;
00311 
00312 
00313         for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00314           if(isamp < sizenew) {
00315             gain_new = (iEB->second)[isamp].gainId();
00316             adc_new = (iEB->second)[isamp].adc();
00317           }
00318           else { adc_new = 0;}
00319 
00320           if(isamp < sizeold) {
00321               gain_old = EB_old[isamp].gainId();
00322               adc_old = EB_old[isamp].adc();
00323           }
00324           else { adc_old = 0;}
00325 
00326           const std::vector<float> pedeStals = GetPedestals(ES,currentID);
00327           const std::vector<float> gainRatios = GetGainRatios(ES,currentID);
00328 
00329           if(adc_new>0 && adc_old>0) {
00330             if(gain_old == gain_new) {  // we're happy - easy case
00331               gain_consensus = gain_old;
00332             }
00333             else {  // lower gain sample has more energy
00334                               
00335 
00336               if(gain_old < gain_new) { // old has higher gain than new, scale to lower gain
00337                 
00338         
00339                 float ratio = gainRatios[gain_new-1]/gainRatios[gain_old-1];
00340                 adc_old = (int) round ((adc_old - pedeStals[gain_old-1]) / ratio + pedeStals[gain_new-1] );  
00341                 gain_consensus = gain_new;
00342               }
00343               else { // scale to old (lower) gain
00344                 float ratio = gainRatios[gain_old-1]/gainRatios[gain_new-1];
00345                 adc_new = (int) round ( (adc_new - pedeStals[gain_new-1]) / ratio+ pedeStals[gain_old-1] );
00346                 gain_consensus = gain_old;
00347               } 
00348             }
00349           }
00350 
00351          
00352           // add values, but don't count pedestals twice
00353           adc_sum = adc_new + adc_old - (int) round (pedeStals[gain_consensus-1]);
00354 
00355 
00356           // if we are now saturating that gain, switch to the next
00357           if (adc_sum> 4096) {
00358             if (gain_consensus<3){
00359 
00360               double ratio = gainRatios[gain_consensus]/gainRatios[gain_consensus-1];
00361               adc_sum = (int) round ((adc_sum - pedeStals[gain_consensus-1])/ ratio + pedeStals[gain_consensus]  )  ;
00362               sw_gain_consensus=++gain_consensus;             
00363             }
00364             else adc_sum = 4096;
00365                 
00366           } 
00367 
00368           // furthermore, make sure we don't decrease our gain once we've switched up
00369           // in case go back 
00370           if (gain_consensus<sw_gain_consensus){
00371 
00372             double ratio = gainRatios[sw_gain_consensus-1]/gainRatios[gain_consensus-1];
00373             adc_sum = (int) round((adc_sum - pedeStals[gain_consensus-1] )/ratio + pedeStals[sw_gain_consensus-1]);
00374             gain_consensus = sw_gain_consensus;
00375           }
00376 
00377           EcalMGPASample sample(adc_sum, gain_consensus);
00378           EB_old.setSample(isamp,sample);  // overwrite old sample, adding new info
00379         } // for sample
00380 
00381 
00382       } // if current = former
00383       else {
00384           if(formerID>0) {
00385             EBdigis->push_back( formerID, EB_old.frame().begin() );
00386           }
00387           //save pointers for next iteration
00388           formerID = currentID;
00389           EB_old = iEB->second;
00390 
00391       }
00392 
00393 
00394       iEBchk = iEB;
00395       if((++iEBchk) == EBDigiStorage_.end()) {  //make sure not to lose the last one
00396             EBdigis->push_back( currentID, (iEB->second).frame().begin()) ;
00397       }
00398     }
00399 
00400     // EE next...
00401 
00402     formerID = 0;
00403     EEDataFrame EE_old;
00404 
00405     EEDigiMap::const_iterator iEEchk;
00406 
00407     for(EEDigiMap::const_iterator iEE  = EEDigiStorage_.begin();
00408         iEE != EEDigiStorage_.end(); iEE++) {
00409 
00410       currentID = iEE->first; 
00411 
00412       if (currentID == formerID) { // we have to add these digis together
00413 
00414         //loop over digi samples in each DataFrame
00415         unsigned int sizenew = (iEE->second).size();
00416         unsigned int sizeold = EE_old.size();
00417 
00418         unsigned int max_samp = std::max(sizenew, sizeold);
00419 
00420         // samples from different events can be of different lengths - sum all
00421         // that overlap.
00422         // check to see if gains match - if not, scale smaller cell down.
00423 
00424         for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00425           if(isamp < sizenew) {
00426             gain_new = (iEE->second)[isamp].gainId();
00427             adc_new = (iEE->second)[isamp].adc();
00428           }
00429           else { adc_new = 0;}
00430 
00431           if(isamp < sizeold) {
00432               gain_old = EE_old[isamp].gainId();
00433               adc_old = EE_old[isamp].adc();
00434           }
00435           else { adc_old = 0;}
00436 
00437           const std::vector<float> pedeStals = GetPedestals(ES,currentID);
00438           const std::vector<float> gainRatios = GetGainRatios(ES,currentID);
00439 
00440           if(adc_new>0 && adc_old>0) {
00441             if(gain_old == gain_new) {  // we're happy - easy case
00442               gain_consensus = gain_old;
00443             }
00444             else {  // lower gain sample has more energy
00445 
00446               if(gain_old < gain_new) { // old has higher gain than new, scale to lower gain
00447                 
00448                 
00449                 float ratio = gainRatios[gain_new-1]/gainRatios[gain_old-1];
00450                 adc_old = (int) round ((adc_old - pedeStals[gain_old-1]) / ratio + pedeStals[gain_new-1] );  
00451                 gain_consensus = gain_new;
00452               }
00453               else { // scale to old (lower) gain
00454                 float ratio = gainRatios[gain_old-1]/gainRatios[gain_new-1];
00455                 adc_new = (int) round ( (adc_new - pedeStals[gain_new-1]) / ratio+ pedeStals[gain_old-1] );
00456                 gain_consensus = gain_old;
00457               } 
00458             }
00459             
00460           }
00461              
00462 
00463           // add values
00464           adc_sum = adc_new + adc_old;
00465           
00466           // if the sum saturates this gain, switch
00467           if (adc_sum> 4096) {
00468             if (gain_consensus<3){
00469               
00470               double ratio = gainRatios[gain_consensus]/gainRatios[gain_consensus-1];
00471               adc_sum = (int) round ((adc_sum - pedeStals[gain_consensus-1])/ ratio + pedeStals[gain_consensus]  )  ;
00472               ++gain_consensus;
00473             }
00474             else adc_sum = 4096;
00475             
00476           } 
00477           
00478           EcalMGPASample sample(adc_sum, gain_consensus);
00479           EE_old.setSample(isamp,sample);
00480         }
00481 
00482       }
00483       else {
00484           if(formerID>0) {
00485             EEdigis->push_back(formerID, EE_old.frame().begin() );
00486              
00487           }
00488           //save pointers for next iteration
00489           formerID = currentID;
00490           EE_old = iEE->second;
00491       }
00492      
00493       iEEchk = iEE;
00494       if((++iEEchk) == EEDigiStorage_.end()) {  //make sure not to lose the last one
00495             EEdigis->push_back(currentID, (iEE->second).frame().begin());
00496       }
00497     }
00498    
00499 
00500     // ES next...
00501 
00502     formerID = 0;
00503     ESDataFrame ES_old;
00504 
00505     ESDigiMap::const_iterator iESchk;
00506 
00507     for(ESDigiMap::const_iterator iES  = ESDigiStorage_.begin();
00508         iES != ESDigiStorage_.end(); iES++) {
00509 
00510       currentID = iES->first; 
00511 
00512       if (currentID == formerID) { // we have to add these digis together
00513 
00514         //loop over digi samples in each DataFrame
00515         unsigned int sizenew = (iES->second).size();
00516         unsigned int sizeold = ES_old.size();
00517         uint16_t rawdat = 0;
00518         unsigned int max_samp = std::max(sizenew, sizeold);
00519 
00520         // samples from different events can be of different lengths - sum all
00521         // that overlap.
00522         // check to see if gains match - if not, scale smaller cell down.
00523 
00524         for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00525           if(isamp < sizenew) {
00526             adc_new = (iES->second)[isamp].adc();
00527             rawdat = (iES->second)[isamp].raw();
00528           }
00529           else { adc_new = 0;}
00530 
00531           if(isamp < sizeold) {
00532               adc_old = ES_old[isamp].adc();
00533               rawdat = ES_old[isamp].raw();
00534           }
00535           else { adc_old = 0;}
00536 
00537           // add values
00538           adc_sum = adc_new + adc_old;
00539           // make data word of gain, rawdata
00540           adc_sum = std::min(adc_sum,4095); //first 12 bits of (uint)
00541           data = adc_sum+(rawdat&0xF000); // data is 14 bit word with gain as MSBs
00542           ES_old.setSample(isamp,data);
00543         }
00544 
00545       }
00546       else {
00547           if(formerID>0) {
00548             ESdigis->push_back(ES_old);
00549           }
00550           //save pointers for next iteration
00551           formerID = currentID;
00552           ES_old = iES->second;
00553       }
00554 
00555       iESchk = iES;
00556       if((++iESchk) == ESDigiStorage_.end()) {  //make sure not to lose the last one
00557         ESdigis->push_back(iES->second);
00558         //      ESDataFrame df( (*ESdigis)->back() );
00559         //for(int isamp=0; isamp<(iES->second).size(); isamp++) {
00560         //  df.setSample(isamp,(iES->second).data[isamp]);
00561         //      }
00562       }
00563     }
00564 
00565 
00566     // done merging
00567 
00568     // put the collection of reconstructed hits in the event   
00569     LogInfo("DataMixingEMDigiWorker") << "total # EB Merged digis: " << EBdigis->size() ;
00570     LogInfo("DataMixingEMDigiWorker") << "total # EE Merged digis: " << EEdigis->size() ;
00571     LogInfo("DataMixingEMDigiWorker") << "total # ES Merged digis: " << ESdigis->size() ;
00572 
00573     e.put( EBdigis, EBDigiCollectionDM_ );
00574     e.put( EEdigis, EEDigiCollectionDM_ );
00575     e.put( ESdigis, ESDigiCollectionDM_ );
00576     
00577     // clear local storage after this event
00578 
00579     EBDigiStorage_.clear();
00580     EEDigiStorage_.clear();
00581     ESDigiStorage_.clear();
00582 
00583   }
00584   const std::vector<float>  DataMixingEMDigiWorker::GetPedestals (const edm::EventSetup& ES, const DetId& detid) {
00585     
00586     std::vector<float> pedeStals(3);
00587 
00588     // get pedestals
00589     edm::ESHandle<EcalPedestals> pedHandle;
00590     ES.get<EcalPedestalsRcd>().get( pedHandle );
00591     
00592     
00593     const EcalPedestalsMap & pedMap = pedHandle.product()->getMap(); // map of pedestals
00594     EcalPedestalsMapIterator pedIter; // pedestal iterator
00595     EcalPedestals::Item aped; // pedestal object for a single xtal
00596     
00597     
00598     pedIter = pedMap.find(detid);
00599     if( pedIter != pedMap.end() ) {
00600       aped = (*pedIter);
00601       pedeStals[0] = aped.mean_x12;
00602       pedeStals[1] = aped.mean_x6;
00603       pedeStals[2] = aped.mean_x1;
00604     } else {
00605       edm::LogError("DataMixingMissingInput") << "Cannot find pedestals";  
00606       pedeStals[0] = 0;
00607       pedeStals[1] = 0;
00608       pedeStals[2] = 0;
00609     }
00610     
00611     
00612     return pedeStals;
00613   }
00614 
00615   const std::vector<float>  DataMixingEMDigiWorker::GetGainRatios(const edm::EventSetup& ES, const DetId& detid) {
00616 
00617     std::vector<float> gainRatios(3);
00618     // get gain ratios  
00619     edm::ESHandle<EcalGainRatios> grHandle;
00620     ES.get<EcalGainRatiosRcd>().get(grHandle);
00621     EcalMGPAGainRatio theRatio= (*grHandle)[detid];
00622     
00623     
00624     gainRatios[0] = 1.;
00625     gainRatios[1] = theRatio.gain12Over6();
00626     gainRatios[2] = theRatio.gain6Over1()  * theRatio.gain12Over6();
00627 
00628     return gainRatios;
00629   }
00630 
00631 
00632 } //edm