CMS 3D CMS Logo

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