CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/CalibCalorimetry/CaloMiscalibTools/plugins/EcalRecHitRecalib.cc

Go to the documentation of this file.
00001 
00002 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
00003 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00004 
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 
00009 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00010 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00011 #include "CalibCalorimetry/CaloMiscalibTools/interface/EcalRecHitRecalib.h"
00012 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 EcalRecHitRecalib::EcalRecHitRecalib(const edm::ParameterSet& iConfig)
00016 {
00017   ecalHitsProducer_ = iConfig.getParameter< std::string > ("ecalRecHitsProducer");
00018   barrelHits_ = iConfig.getParameter< std::string > ("barrelHitCollection");
00019   endcapHits_ = iConfig.getParameter< std::string > ("endcapHitCollection");
00020   RecalibBarrelHits_ = iConfig.getParameter< std::string > ("RecalibBarrelHitCollection");
00021   RecalibEndcapHits_ = iConfig.getParameter< std::string > ("RecalibEndcapHitCollection");
00022   refactor_ = iConfig.getUntrackedParameter<double> ("Refactor",(double)1);
00023   refactor_mean_ = iConfig.getUntrackedParameter<double> ("Refactor_mean",(double)1);
00024 
00025   //register your products
00026   produces< EBRecHitCollection >(RecalibBarrelHits_);
00027   produces< EERecHitCollection >(RecalibEndcapHits_);
00028 }
00029 
00030 
00031 EcalRecHitRecalib::~EcalRecHitRecalib()
00032 {
00033  
00034 
00035 }
00036 
00037 
00038 // ------------ method called to produce the data  ------------
00039 void
00040 EcalRecHitRecalib::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00041 {
00042   using namespace edm;
00043   using namespace std;
00044 
00045   Handle<EBRecHitCollection> barrelRecHitsHandle;
00046   Handle<EERecHitCollection> endcapRecHitsHandle;
00047 
00048   const EBRecHitCollection*  EBRecHits = 0;
00049   const EERecHitCollection*  EERecHits = 0; 
00050  
00051   iEvent.getByLabel(ecalHitsProducer_,barrelHits_,barrelRecHitsHandle);
00052   if (!barrelRecHitsHandle.isValid()) {
00053     LogDebug("") << "EcalREcHitMiscalib: Error! can't get product!" << std::endl;
00054   } else {
00055     EBRecHits = barrelRecHitsHandle.product(); // get a ptr to the product
00056   }
00057 
00058   iEvent.getByLabel(ecalHitsProducer_,endcapHits_,endcapRecHitsHandle);
00059   if (!endcapRecHitsHandle.isValid()) {
00060     LogDebug("") << "EcalREcHitMiscalib: Error! can't get product!" << std::endl;
00061   } else {
00062     EERecHits = endcapRecHitsHandle.product(); // get a ptr to the product
00063   }
00064 
00065   //Create empty output collections
00066   std::auto_ptr< EBRecHitCollection > RecalibEBRecHitCollection( new EBRecHitCollection );
00067   std::auto_ptr< EERecHitCollection > RecalibEERecHitCollection( new EERecHitCollection );
00068 
00069 
00070   // Intercalib constants
00071   edm::ESHandle<EcalIntercalibConstants> pIcal;
00072   iSetup.get<EcalIntercalibConstantsRcd>().get(pIcal);
00073   const EcalIntercalibConstants* ical = pIcal.product();
00074 
00075   if(EBRecHits)
00076     {
00077 
00078        //loop on all EcalRecHits (barrel)
00079       EBRecHitCollection::const_iterator itb;
00080       for (itb=EBRecHits->begin(); itb!=EBRecHits->end(); itb++) {
00081         
00082         // find intercalib constant for this xtal
00083         EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
00084         EcalIntercalibConstant icalconst = -1;
00085 
00086         if( icalit!=ical->getMap().end() ){
00087           icalconst = (*icalit);
00088           // edm::LogDebug("EcalRecHitRecalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;
00089 
00090         } else {
00091           edm::LogError("EcalRecHitRecalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00092               ;
00093           }
00094           
00095           // make the rechit with rescaled energy and put in the output collection
00096         icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)
00097         EcalRecHit aHit(itb->id(),itb->energy()*icalconst,itb->time());
00098           
00099         RecalibEBRecHitCollection->push_back( aHit);
00100       }
00101     }
00102 
00103   if(EERecHits)
00104     {
00105 
00106        //loop on all EcalRecHits (barrel)
00107       EERecHitCollection::const_iterator ite;
00108       for (ite=EERecHits->begin(); ite!=EERecHits->end(); ite++) {
00109         
00110         // find intercalib constant for this xtal
00111         EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(ite->id().rawId());
00112         EcalIntercalibConstant icalconst = -1;
00113 
00114         if( icalit!=ical->getMap().end() ){
00115           icalconst = (*icalit);
00116           // edm:: LogDebug("EcalRecHitRecalib") << "Found intercalib for xtal " << EEDetId(ite->id()) << " " << icalconst ;
00117           } else {
00118             edm::LogError("EcalRecHitRecalib") << "No intercalib const found for xtal " << EEDetId(ite->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00119               ;
00120           }
00121           
00122           // make the rechit with rescaled energy and put in the output collection
00123         
00124         icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)
00125         EcalRecHit aHit(ite->id(),ite->energy()*icalconst,ite->time());
00126         
00127         RecalibEERecHitCollection->push_back( aHit);
00128       }
00129     }
00130 
00131 
00132   //Put Recalibrated rechit in the event
00133   iEvent.put( RecalibEBRecHitCollection, RecalibBarrelHits_);
00134   iEvent.put( RecalibEERecHitCollection, RecalibEndcapHits_);
00135   
00136 }
00137