CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/CalibCalorimetry/CaloMiscalibTools/plugins/HcalRecHitRecalib.cc

Go to the documentation of this file.
00001 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00002 #include "CalibCalorimetry/CaloMiscalibTools/interface/HcalRecHitRecalib.h"
00003 
00004 #include "DataFormats/Common/interface/Handle.h"
00005 #include "FWCore/ParameterSet/interface/FileInPath.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 
00008 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "CalibCalorimetry/CaloMiscalibTools/interface/MiscalibReaderFromXMLHcal.h"
00011 
00012 HcalRecHitRecalib::HcalRecHitRecalib(const edm::ParameterSet& iConfig)
00013 {
00014 
00015   hbheLabel_ = iConfig.getParameter<edm::InputTag>("hbheInput");
00016   hoLabel_ = iConfig.getParameter<edm::InputTag>("hoInput");
00017   hfLabel_ = iConfig.getParameter<edm::InputTag>("hfInput");
00018 
00019 //   HBHEHitsProducer_ = iConfig.getParameter< std::string > ("HBHERecHitsProducer");
00020 //   HOHitsProducer_ = iConfig.getParameter< std::string > ("HERecHitsProducer");
00021 //   HFHitsProducer_ = iConfig.getParameter< std::string > ("HERecHitsProducer");
00022 //   HBHEHits_ = iConfig.getParameter< std::string > ("HBHEHitCollection");
00023 //   HFHits_ = iConfig.getParameter< std::string > ("HFHitCollection");
00024 //   HOHits_ = iConfig.getParameter< std::string > ("HOHitCollection");
00025 
00026   RecalibHBHEHits_ = iConfig.getParameter< std::string > ("RecalibHBHEHitCollection");
00027   RecalibHFHits_ = iConfig.getParameter< std::string > ("RecalibHFHitCollection");
00028   RecalibHOHits_ = iConfig.getParameter< std::string > ("RecalibHOHitCollection");
00029 
00030   refactor_ = iConfig.getUntrackedParameter<double> ("Refactor",(double)1);
00031   refactor_mean_ = iConfig.getUntrackedParameter<double> ("Refactor_mean",(double)1);
00032 
00033   //register your products
00034   produces< HBHERecHitCollection >(RecalibHBHEHits_);
00035   produces< HFRecHitCollection >(RecalibHFHits_);
00036   produces< HORecHitCollection >(RecalibHOHits_);
00037 
00038   // here read them from xml (particular to HCAL)
00039 
00040   hcalfileinpath_=iConfig.getUntrackedParameter<std::string> ("fileNameHcal","");
00041   edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);
00042 
00043   hcalfile_=hcalfiletmp.fullPath();
00044 }
00045 
00046 
00047 HcalRecHitRecalib::~HcalRecHitRecalib()
00048 {
00049  
00050 
00051 }
00052 
00053 void
00054 HcalRecHitRecalib::beginRun(const edm::Run&, const edm::EventSetup& iSetup)
00055 {
00056   edm::ESHandle<HcalTopology> topology;
00057   iSetup.get<IdealGeometryRecord>().get( topology );
00058   
00059   mapHcal_.prefillMap(*topology);
00060 
00061   MiscalibReaderFromXMLHcal hcalreader_(mapHcal_);
00062   if(!hcalfile_.empty()) hcalreader_.parseXMLMiscalibFile(hcalfile_);
00063   mapHcal_.print();
00064 }
00065 
00066 // ------------ method called to produce the data  ------------
00067 void
00068 HcalRecHitRecalib::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00069 {
00070   using namespace edm;
00071   using namespace std;
00072 
00073   Handle<HBHERecHitCollection> HBHERecHitsHandle;
00074   Handle<HFRecHitCollection> HFRecHitsHandle;
00075   Handle<HORecHitCollection> HORecHitsHandle;
00076 
00077   const HBHERecHitCollection*  HBHERecHits = 0;
00078   const HFRecHitCollection*  HFRecHits = 0;
00079   const HORecHitCollection*  HORecHits = 0;
00080 
00081   iEvent.getByLabel(hbheLabel_,HBHERecHitsHandle);
00082   if (!HBHERecHitsHandle.isValid()) {
00083     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
00084   } else {
00085     HBHERecHits = HBHERecHitsHandle.product(); // get a ptr to the product
00086   }
00087 
00088   iEvent.getByLabel(hoLabel_,HORecHitsHandle);
00089   if (!HORecHitsHandle.isValid()) {
00090     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
00091   } else {
00092     HORecHits = HORecHitsHandle.product(); // get a ptr to the product
00093   }
00094 
00095   iEvent.getByLabel(hfLabel_,HFRecHitsHandle);
00096   if (!HFRecHitsHandle.isValid()) {
00097     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
00098   } else {
00099     HFRecHits = HFRecHitsHandle.product(); // get a ptr to the product
00100   }
00101 
00102 
00103 //     iEvent.getByLabel(HBHEHitsProducer_,HBHEHits_,HBHERecHitsHandle);
00104 //     HBHERecHits = HBHERecHitsHandle.product(); // get a ptr to the product
00105 
00106 //     iEvent.getByLabel(HFHitsProducer_,HFHits_,HFRecHitsHandle);
00107 //     HFRecHits = HFRecHitsHandle.product(); // get a ptr to the product
00108 
00109 //     iEvent.getByLabel(HOHitsProducer_,HOHits_,HORecHitsHandle);
00110 //     HORecHits = HORecHitsHandle.product(); // get a ptr to the product
00111 
00112 
00113   //Create empty output collections
00114   std::auto_ptr< HBHERecHitCollection > RecalibHBHERecHitCollection( new HBHERecHitCollection );
00115   std::auto_ptr< HFRecHitCollection > RecalibHFRecHitCollection( new HFRecHitCollection );
00116   std::auto_ptr< HORecHitCollection > RecalibHORecHitCollection( new HORecHitCollection );
00117 
00118   // Intercalib constants
00119   //  edm::ESHandle<EcalIntercalibConstants> pIcal;
00120   //  iSetup.get<EcalIntercalibConstantsRcd>().get(pIcal);
00121   //  const EcalIntercalibConstants* ical = pIcal.product();
00122 
00123   if(HBHERecHits)
00124     {
00125 
00126        //loop on all EcalRecHits (barrel)
00127       HBHERecHitCollection::const_iterator itHBHE;
00128       for (itHBHE=HBHERecHits->begin(); itHBHE!=HBHERecHits->end(); itHBHE++) {
00129         
00130         // find intercalib constant for this cell
00131 
00132         //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
00133         //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;
00134 
00135         //      if( icalit!=ical->getMap().end() ){
00136         //        icalconst = icalit->second;
00137           // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;
00138 
00139         //      } else {
00140         //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00141         //              ;
00142         //          }
00143           
00144         float icalconst=(mapHcal_.get().find(itHBHE->id().rawId()))->second;
00145           // make the rechit with rescaled energy and put in the output collection
00146 
00147         icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
00148         HBHERecHit aHit(itHBHE->id(),itHBHE->energy()*icalconst,itHBHE->time());
00149         
00150         RecalibHBHERecHitCollection->push_back( aHit);
00151       }
00152     }
00153 
00154   if(HFRecHits)
00155     {
00156 
00157        //loop on all EcalRecHits (barrel)
00158       HFRecHitCollection::const_iterator itHF;
00159       for (itHF=HFRecHits->begin(); itHF!=HFRecHits->end(); itHF++) {
00160         
00161         // find intercalib constant for this cell
00162 
00163         //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
00164         //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;
00165 
00166         //      if( icalit!=ical->getMap().end() ){
00167         //        icalconst = icalit->second;
00168           // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;
00169 
00170         //      } else {
00171         //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00172         //              ;
00173         //          }
00174           
00175           // make the rechit with rescaled energy and put in the output collection
00176         
00177         float icalconst=(mapHcal_.get().find(itHF->id().rawId()))->second;
00178         icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
00179         HFRecHit aHit(itHF->id(),itHF->energy()*icalconst,itHF->time());
00180         
00181         RecalibHFRecHitCollection->push_back( aHit);
00182       }
00183     }
00184 
00185   if(HORecHits)
00186     {
00187 
00188        //loop on all EcalRecHits (barrel)
00189       HORecHitCollection::const_iterator itHO;
00190       for (itHO=HORecHits->begin(); itHO!=HORecHits->end(); itHO++) {
00191         
00192         // find intercalib constant for this cell
00193 
00194         //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
00195         //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;
00196 
00197         //      if( icalit!=ical->getMap().end() ){
00198         //        icalconst = icalit->second;
00199           // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;
00200 
00201         //      } else {
00202         //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
00203         //              ;
00204         //          }
00205           
00206           // make the rechit with rescaled energy and put in the output collection
00207 
00208         float icalconst=(mapHcal_.get().find(itHO->id().rawId()))->second;
00209         icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
00210         HORecHit aHit(itHO->id(),itHO->energy()*icalconst,itHO->time());
00211           
00212           RecalibHORecHitCollection->push_back( aHit);
00213       }
00214     }
00215 
00216 
00217   //Put Recalibrated rechit in the event
00218   iEvent.put( RecalibHBHERecHitCollection, RecalibHBHEHits_);
00219   iEvent.put( RecalibHFRecHitCollection, RecalibHFHits_);
00220   iEvent.put( RecalibHORecHitCollection, RecalibHOHits_);
00221 }
00222