CMS 3D CMS Logo

EcalRecalibRecHitProducer.cc

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

Generated on Tue Jun 9 17:43:47 2009 for CMSSW by  doxygen 1.5.4