CMS 3D CMS Logo

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