CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Calibration/HcalCalibAlgos/src/HitReCalibrator.cc

Go to the documentation of this file.
00001 #include "Calibration/HcalCalibAlgos/src/HitReCalibrator.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00004 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00005 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00006 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00007 #include "RecoTracker/TrackProducer/interface/TrackProducerBase.h"
00008 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00009 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00010 #include "CondFormats/HcalObjects/interface/HcalRespCorrs.h"
00011 #include "CondFormats/DataRecord/interface/HcalRespCorrsRcd.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 
00014 
00015 using namespace edm;
00016 using namespace std;
00017 using namespace reco;
00018 
00019 namespace cms
00020 {
00021 
00022 HitReCalibrator::HitReCalibrator(const edm::ParameterSet& iConfig)
00023 {
00024    hbheInput_ = iConfig.getParameter<edm::InputTag>("hbheInput");
00025    hoInput_ = iConfig.getParameter<edm::InputTag>("hoInput");
00026    hfInput_ = iConfig.getParameter<edm::InputTag>("hfInput"); 
00027    allowMissingInputs_ = true;
00028 //register your products
00029 
00030    produces<HBHERecHitCollection>("DiJetsHBHEReRecHitCollection");
00031    produces<HORecHitCollection>("DiJetsHOReRecHitCollection");
00032    produces<HFRecHitCollection>("DiJetsHFReRecHitCollection");
00033 
00034 }
00035 void HitReCalibrator::beginJob()
00036 {
00037 }
00038 
00039 HitReCalibrator::~HitReCalibrator()
00040 {
00041 
00042 }
00043 
00044 
00045 // ------------ method called to produce the data  ------------
00046 void
00047 HitReCalibrator::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00048 {
00049 
00050   std::auto_ptr<HBHERecHitCollection> miniDiJetsHBHERecHitCollection(new HBHERecHitCollection);
00051   std::auto_ptr<HORecHitCollection> miniDiJetsHORecHitCollection(new HORecHitCollection);
00052   std::auto_ptr<HFRecHitCollection> miniDiJetsHFRecHitCollection(new HFRecHitCollection);
00053 
00054  
00055   edm::ESHandle <HcalRespCorrs> recalibCorrs;
00056   iSetup.get<HcalRespCorrsRcd>().get("recalibrate",recalibCorrs);
00057   const HcalRespCorrs* jetRecalib = recalibCorrs.product();
00058  
00059    
00060   try {
00061   edm::Handle<HBHERecHitCollection> hbhe;
00062   iEvent.getByLabel(hbheInput_,hbhe);
00063   const HBHERecHitCollection Hithbhe = *(hbhe.product());
00064   for(HBHERecHitCollection::const_iterator hbheItr=Hithbhe.begin(); hbheItr!=Hithbhe.end(); hbheItr++)
00065         {
00066            DetId id = hbheItr->detid();
00067            float recal; 
00068            if (jetRecalib->exists(id))
00069              recal = jetRecalib->getValues(id)->getValue();
00070            else recal = 1.; 
00071            float energy = hbheItr->energy(); 
00072            float time = hbheItr->time();
00073            HBHERecHit* hbhehit = new HBHERecHit(id,recal*energy,time);  
00074            miniDiJetsHBHERecHitCollection->push_back(*hbhehit);
00075         }
00076   } catch (cms::Exception& e) { // can't find it!
00077   if (!allowMissingInputs_) {cout<<"No HBHE collection "<<endl; throw e;}
00078   }
00079 
00080   try{  
00081   edm::Handle<HORecHitCollection> ho;
00082   iEvent.getByLabel(hoInput_,ho);
00083   const HORecHitCollection Hitho = *(ho.product());
00084   for(HORecHitCollection::const_iterator hoItr=Hitho.begin(); hoItr!=Hitho.end(); hoItr++)
00085         {
00086            DetId id = hoItr->detid();
00087            float recal; 
00088            if (jetRecalib->exists(id))
00089              recal = jetRecalib->getValues(id)->getValue();
00090            else recal = 1.; 
00091            float energy = hoItr->energy(); 
00092            float time = hoItr->time();
00093            HORecHit* hohit = new HORecHit(id,recal*energy,time); 
00094            miniDiJetsHORecHitCollection->push_back(*hohit);
00095         }
00096   } catch (cms::Exception& e) { // can't find it!
00097   if (!allowMissingInputs_) {cout<<" No HO collection "<<endl; throw e;}
00098   }
00099 
00100   try {
00101   edm::Handle<HFRecHitCollection> hf;
00102   iEvent.getByLabel(hfInput_,hf);
00103   const HFRecHitCollection Hithf = *(hf.product());
00104   for(HFRecHitCollection::const_iterator hfItr=Hithf.begin(); hfItr!=Hithf.end(); hfItr++)
00105       {
00106            DetId id = hfItr->detid();
00107            float recal; 
00108            if (jetRecalib->exists(id))
00109              recal = jetRecalib->getValues(id)->getValue();
00110            else recal = 1.; 
00111            float energy = hfItr->energy(); 
00112            float time = hfItr->time();
00113            HFRecHit* hfhit = new HFRecHit(id,recal*energy,time); 
00114            miniDiJetsHFRecHitCollection->push_back(*hfhit);
00115       }
00116   } catch (cms::Exception& e) { // can't find it!
00117   if (!allowMissingInputs_) throw e;
00118   }
00119 
00120   //Put selected information in the event
00121 
00122   iEvent.put( miniDiJetsHBHERecHitCollection, "DiJetsHBHEReRecHitCollection");
00123 
00124   iEvent.put( miniDiJetsHORecHitCollection, "DiJetsHOReRecHitCollection");
00125 
00126   iEvent.put( miniDiJetsHFRecHitCollection, "DiJetsHFReRecHitCollection");
00127 
00128 }
00129 }