CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/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         if ( EBRecHitCollection_.label() != "" ) {
00074                 evt.getByLabel( EBRecHitCollection_, pEBRecHits);
00075                 if ( pEBRecHits.isValid() ) {
00076                         EBRecHits = pEBRecHits.product(); // get a ptr to the product
00077 #ifdef DEBUG
00078                         LogDebug("EcalRecHitDebug") << "total # EB rechits to be re-calibrated: " << EBRecHits->size();
00079 #endif
00080                 } else {
00081                         edm::LogError("EcalRecHitError") << "Error! can't get the product " << EBRecHitCollection_.label() ;
00082                 }
00083         }
00084 
00085         //        if ( EERecHitCollection_.label() != "" && EERecHitCollection_.instance() != "" ) {
00086         if ( EERecHitCollection_.label() != ""  ) {
00087                 evt.getByLabel( EERecHitCollection_, pEERecHits);
00088                 if ( pEERecHits.isValid() ) {
00089                         EERecHits = pEERecHits.product(); // get a ptr to the product
00090 #ifdef DEBUG
00091                         LogDebug("EcalRecHitDebug") << "total # EE uncalibrated rechits to be re-calibrated: " << EERecHits->size();
00092 #endif
00093                 } else {
00094                         edm::LogError("EcalRecHitError") << "Error! can't get the product " << EERecHitCollection_.label() ;
00095                 }
00096         }
00097 
00098         // collection of rechits to put in the event
00099         std::auto_ptr< EBRecHitCollection > EBRecalibRecHits( new EBRecHitCollection );
00100         std::auto_ptr< EERecHitCollection > EERecalibRecHits( new EERecHitCollection );
00101 
00102         // now fetch all conditions we need to make rechits
00103         // ADC to GeV constant
00104         edm::ESHandle<EcalADCToGeVConstant> pAgc;
00105         const EcalADCToGeVConstant *agc = 0;
00106         float agc_eb = 1.;
00107         float agc_ee = 1.;
00108         if (doEnergyScale_) {
00109                 es.get<EcalADCToGeVConstantRcd>().get(pAgc);
00110                 agc = pAgc.product();
00111                 // use this value in the algorithm
00112                 agc_eb = float(agc->getEBValue());
00113                 agc_ee = float(agc->getEEValue());
00114         }
00115         // Intercalib constants
00116         edm::ESHandle<EcalIntercalibConstants> pIcal;
00117         const EcalIntercalibConstants *ical = 0;
00118         if (doIntercalib_) {
00119                 es.get<EcalIntercalibConstantsRcd>().get(pIcal);
00120                 ical = pIcal.product();
00121         }
00122         // Laser corrections
00123         edm::ESHandle<EcalLaserDbService> pLaser;
00124         es.get<EcalLaserDbRecord>().get( pLaser );
00125 
00126         if (EBRecHits) {
00127                 // loop over uncalibrated rechits to make calibrated ones
00128                 for(EBRecHitCollection::const_iterator it  = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
00129 
00130                         EcalIntercalibConstant icalconst = 1.;
00131                         if (doIntercalib_) {
00132                                 // find intercalib constant for this xtal
00133                                 const EcalIntercalibConstantMap &icalMap = ical->getMap();
00134                                 EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
00135                                 if( icalit!=icalMap.end() ){
00136                                         icalconst = (*icalit);
00137                                 } else {
00138                                         edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00139                                                 ;
00140                                 }
00141                         }
00142                         // get laser coefficient
00143                         float lasercalib = 1;
00144                         if (doLaserCorrections_) {
00145                                 lasercalib = pLaser->getLaserCorrection( EBDetId(it->id()), evt.time() );
00146                         }
00147 
00148                         // make the rechit and put in the output collection
00149                         // must implement op= for EcalRecHit
00150                         EcalRecHit aHit( (*it).id(), (*it).energy() * agc_eb * icalconst * lasercalib, (*it).time() );
00151                         EBRecalibRecHits->push_back( aHit );
00152                 }
00153         }
00154 
00155         if (EERecHits)
00156         {
00157                 // loop over uncalibrated rechits to make calibrated ones
00158                 for(EERecHitCollection::const_iterator it  = EERecHits->begin();
00159                                 it != EERecHits->end(); ++it) {
00160 
00161                         // find intercalib constant for this xtal
00162                         EcalIntercalibConstant icalconst = 1.;
00163                         if (doIntercalib_) {
00164                                 const EcalIntercalibConstantMap &icalMap = ical->getMap();
00165                                 EcalIntercalibConstantMap::const_iterator icalit=icalMap.find(it->id());
00166                                 if( icalit!=icalMap.end() ) {
00167                                         icalconst = (*icalit);
00168                                 } else {
00169                                         edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EEDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? ";
00170                                 }
00171                         }
00172                         // get laser coefficient
00173                         float lasercalib = 1;
00174                         if (doLaserCorrections_) {
00175                                 lasercalib = pLaser->getLaserCorrection( EEDetId(it->id()), evt.time() );
00176                         }
00177 
00178                         // make the rechit and put in the output collection
00179                         // must implement op= for EcalRecHit
00180                         //EcalRecHit aHit( EEalgo_->makeRecHit(*it, icalconst * lasercalib) );
00181                         EcalRecHit aHit( (*it).id(), (*it).energy() * agc_ee * icalconst * lasercalib, (*it).time() );
00182                         EERecalibRecHits->push_back( aHit );
00183                 }
00184         }
00185         // put the collection of recunstructed hits in the event   
00186         LogInfo("EcalRecalibRecHitInfo") << "total # EB re-calibrated rechits: " << EBRecalibRecHits->size();
00187         LogInfo("EcalRecalibRecHitInfo") << "total # EE re-calibrated rechits: " << EERecalibRecHits->size();
00188 
00189         evt.put( EBRecalibRecHits, EBRecalibRecHitCollection_ );
00190         evt.put( EERecalibRecHits, EERecalibRecHitCollection_ );
00191 }
00192 
00193 #include "FWCore/Framework/interface/MakerMacros.h"
00194 DEFINE_FWK_MODULE( EcalRecalibRecHitProducer );