CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/CalibCalorimetry/CaloMiscalibTools/plugins/HcalRecHitRecalib.cc

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