CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalCalo/EcalRecProducers/plugins/EcalRecalibRecHitProducer.cc

Go to the documentation of this file.
00001 
00011 #include "RecoLocalCalo/EcalRecProducers/plugins/EcalRecalibRecHitProducer.h"
00012 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalRecHitSimpleAlgo.h"
00013 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
00014 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00015 #include "CondFormats/EcalObjects/interface/EcalADCToGeVConstant.h"
00016 #include "CondFormats/DataRecord/interface/EcalADCToGeVConstantRcd.h"
00017 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbService.h"
00018 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbRecord.h"
00019 
00020 #include "DataFormats/Common/interface/Handle.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 
00025 #include <iostream>
00026 #include <cmath>
00027 #include <vector>
00028 
00029 //#include "CondFormats/EcalObjects/interface/EcalPedestals.h"
00030 //#include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
00031 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00032 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
00033 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00034 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00035 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00036 
00037 
00038 
00039 EcalRecalibRecHitProducer::EcalRecalibRecHitProducer(const edm::ParameterSet& ps) {
00040 
00041    EBRecHitCollection_ = ps.getParameter<edm::InputTag>("EBRecHitCollection");
00042    EERecHitCollection_ = ps.getParameter<edm::InputTag>("EERecHitCollection");
00043    EBRecalibRecHitCollection_        = ps.getParameter<std::string>("EBRecalibRecHitCollection");
00044    EERecalibRecHitCollection_        = ps.getParameter<std::string>("EERecalibRecHitCollection");
00045    doEnergyScale_             = ps.getParameter<bool>("doEnergyScale");
00046    doIntercalib_              = ps.getParameter<bool>("doIntercalib");
00047    doLaserCorrections_        = ps.getParameter<bool>("doLaserCorrections");
00048 
00049    EBalgo_ = new EcalRecHitSimpleAlgo();
00050    EEalgo_ = new EcalRecHitSimpleAlgo();
00051 
00052    produces< EBRecHitCollection >(EBRecalibRecHitCollection_);
00053    produces< EERecHitCollection >(EERecalibRecHitCollection_);
00054 }
00055 
00056 EcalRecalibRecHitProducer::~EcalRecalibRecHitProducer() {
00057 
00058   if (EBalgo_) delete EBalgo_;
00059   if (EEalgo_) delete EEalgo_;
00060 
00061 }
00062 
00063 void EcalRecalibRecHitProducer::produce(edm::Event& evt, const edm::EventSetup& es)
00064 {
00065         using namespace edm;
00066         Handle< EBRecHitCollection > pEBRecHits;
00067         Handle< EERecHitCollection > pEERecHits;
00068 
00069         const EBRecHitCollection*  EBRecHits = 0;
00070         const EERecHitCollection*  EERecHits = 0; 
00071 
00072         if ( EBRecHitCollection_.label() != "" && EBRecHitCollection_.instance() != "" ) {
00073                 evt.getByLabel( EBRecHitCollection_, pEBRecHits);
00074                 if ( pEBRecHits.isValid() ) {
00075                         EBRecHits = pEBRecHits.product(); // get a ptr to the product
00076 #ifdef DEBUG
00077                         LogDebug("EcalRecHitDebug") << "total # EB rechits to be re-calibrated: " << EBRecHits->size();
00078 #endif
00079                 } else {
00080                         edm::LogError("EcalRecHitError") << "Error! can't get the product " << EBRecHitCollection_.label() ;
00081                 }
00082         }
00083 
00084         if ( EERecHitCollection_.label() != "" && EERecHitCollection_.instance() != "" ) {
00085                 evt.getByLabel( EERecHitCollection_, pEERecHits);
00086                 if ( pEERecHits.isValid() ) {
00087                         EERecHits = pEERecHits.product(); // get a ptr to the product
00088 #ifdef DEBUG
00089                         LogDebug("EcalRecHitDebug") << "total # EE uncalibrated rechits to be re-calibrated: " << EERecHits->size();
00090 #endif
00091                 } else {
00092                         edm::LogError("EcalRecHitError") << "Error! can't get the product " << EERecHitCollection_.label() ;
00093                 }
00094         }
00095 
00096         // collection of rechits to put in the event
00097         std::auto_ptr< EBRecHitCollection > EBRecalibRecHits( new EBRecHitCollection );
00098         std::auto_ptr< EERecHitCollection > EERecalibRecHits( new EERecHitCollection );
00099 
00100         // now fetch all conditions we need to make rechits
00101         // ADC to GeV constant
00102         edm::ESHandle<EcalADCToGeVConstant> pAgc;
00103         const EcalADCToGeVConstant *agc = 0;
00104         float agc_eb = 1.;
00105         float agc_ee = 1.;
00106         if (doEnergyScale_) {
00107                 es.get<EcalADCToGeVConstantRcd>().get(pAgc);
00108                 agc = pAgc.product();
00109                 // use this value in the algorithm
00110                 agc_eb = float(agc->getEBValue());
00111                 agc_ee = float(agc->getEEValue());
00112         }
00113         // Intercalib constants
00114         edm::ESHandle<EcalIntercalibConstants> pIcal;
00115         const EcalIntercalibConstants *ical = 0;
00116         if (doIntercalib_) {
00117                 es.get<EcalIntercalibConstantsRcd>().get(pIcal);
00118                 ical = pIcal.product();
00119         }
00120         // Laser corrections
00121         edm::ESHandle<EcalLaserDbService> pLaser;
00122         es.get<EcalLaserDbRecord>().get( pLaser );
00123 
00124         if (EBRecHits) {
00125                 // loop over uncalibrated rechits to make calibrated ones
00126                 for(EBRecHitCollection::const_iterator it  = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
00127 
00128                         EcalIntercalibConstant icalconst = 1.;
00129                         if (doIntercalib_) {
00130                                 // find intercalib constant for this xtal
00131                                 const EcalIntercalibConstantMap &icalMap = ical->getMap();
00132                                 EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
00133                                 if( icalit!=icalMap.end() ){
00134                                         icalconst = (*icalit);
00135                                 } else {
00136                                         edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00137                                                 ;
00138                                 }
00139                         }
00140                         // get laser coefficient
00141                         float lasercalib = 1;
00142                         if (doLaserCorrections_) {
00143                                 lasercalib = pLaser->getLaserCorrection( EBDetId(it->id()), evt.time() );
00144                         }
00145 
00146                         // make the rechit and put in the output collection
00147                         // must implement op= for EcalRecHit
00148                         EcalRecHit aHit( (*it).id(), (*it).energy() * agc_eb * icalconst * lasercalib, (*it).time() );
00149                         EBRecalibRecHits->push_back( aHit );
00150                 }
00151         }
00152 
00153         if (EERecHits)
00154         {
00155                 // loop over uncalibrated rechits to make calibrated ones
00156                 for(EERecHitCollection::const_iterator it  = EERecHits->begin();
00157                                 it != EERecHits->end(); ++it) {
00158 
00159                         // find intercalib constant for this xtal
00160                         EcalIntercalibConstant icalconst = 1.;
00161                         if (doIntercalib_) {
00162                                 const EcalIntercalibConstantMap &icalMap = ical->getMap();
00163                                 EcalIntercalibConstantMap::const_iterator icalit=icalMap.find(it->id());
00164                                 if( icalit!=icalMap.end() ) {
00165                                         icalconst = (*icalit);
00166                                 } else {
00167                                         edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EEDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? ";
00168                                 }
00169                         }
00170                         // get laser coefficient
00171                         float lasercalib = 1;
00172                         if (doLaserCorrections_) {
00173                                 lasercalib = pLaser->getLaserCorrection( EEDetId(it->id()), evt.time() );
00174                         }
00175 
00176                         // make the rechit and put in the output collection
00177                         // must implement op= for EcalRecHit
00178                         //EcalRecHit aHit( EEalgo_->makeRecHit(*it, icalconst * lasercalib) );
00179                         EcalRecHit aHit( (*it).id(), (*it).energy() * agc_ee * icalconst * lasercalib, (*it).time() );
00180                         EERecalibRecHits->push_back( aHit );
00181                 }
00182         }
00183         // put the collection of recunstructed hits in the event   
00184         LogInfo("EcalRecalibRecHitInfo") << "total # EB re-calibrated rechits: " << EBRecalibRecHits->size();
00185         LogInfo("EcalRecalibRecHitInfo") << "total # EE re-calibrated rechits: " << EERecalibRecHits->size();
00186 
00187         evt.put( EBRecalibRecHits, EBRecalibRecHitCollection_ );
00188         evt.put( EERecalibRecHits, EERecalibRecHitCollection_ );
00189 }
00190 
00191 #include "FWCore/Framework/interface/MakerMacros.h"
00192 DEFINE_FWK_MODULE( EcalRecalibRecHitProducer );