Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:33

0001 #include "CalibCalorimetry/CaloMiscalibTools/interface/HcalRecHitRecalib.h"
0002 
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0005 #include "FWCore/ParameterSet/interface/FileInPath.h"
0006 
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "CalibCalorimetry/CaloMiscalibTools/interface/MiscalibReaderFromXMLHcal.h"
0009 
0010 HcalRecHitRecalib::HcalRecHitRecalib(const edm::ParameterSet& iConfig)
0011     : tok_hbhe_(consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"))),
0012       tok_ho_(consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"))),
0013       tok_hf_(consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"))),
0014       topologyToken_(esConsumes<edm::Transition::BeginRun>()),
0015       recalibHBHEHits_(iConfig.getParameter<std::string>("RecalibHBHEHitCollection")),
0016       recalibHFHits_(iConfig.getParameter<std::string>("RecalibHFHitCollection")),
0017       recalibHOHits_(iConfig.getParameter<std::string>("RecalibHOHitCollection")),
0018       hcalfileinpath_(iConfig.getUntrackedParameter<std::string>("fileNameHcal", "")),
0019       refactor_(iConfig.getUntrackedParameter<double>("Refactor", (double)1)),
0020       refactor_mean_(iConfig.getUntrackedParameter<double>("Refactor_mean", (double)1)) {
0021   //register your products
0022   produces<HBHERecHitCollection>(recalibHBHEHits_);
0023   produces<HFRecHitCollection>(recalibHFHits_);
0024   produces<HORecHitCollection>(recalibHOHits_);
0025 
0026   // here read them from xml (particular to HCAL)
0027   edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/" + hcalfileinpath_);
0028   hcalfile_ = hcalfiletmp.fullPath();
0029 }
0030 
0031 HcalRecHitRecalib::~HcalRecHitRecalib() {}
0032 
0033 void HcalRecHitRecalib::beginRun(const edm::Run&, const edm::EventSetup& iSetup) {
0034   const HcalTopology& topology = iSetup.getData(topologyToken_);
0035 
0036   mapHcal_.prefillMap(topology);
0037 
0038   MiscalibReaderFromXMLHcal hcalreader_(mapHcal_);
0039   if (!hcalfile_.empty())
0040     hcalreader_.parseXMLMiscalibFile(hcalfile_);
0041   mapHcal_.print();
0042 }
0043 
0044 // ------------ method called to produce the data  ------------
0045 void HcalRecHitRecalib::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0046   using namespace edm;
0047   using namespace std;
0048 
0049   Handle<HBHERecHitCollection> HBHERecHitsHandle;
0050   Handle<HFRecHitCollection> HFRecHitsHandle;
0051   Handle<HORecHitCollection> HORecHitsHandle;
0052 
0053   const HBHERecHitCollection* HBHERecHits = nullptr;
0054   const HFRecHitCollection* HFRecHits = nullptr;
0055   const HORecHitCollection* HORecHits = nullptr;
0056 
0057   iEvent.getByToken(tok_hbhe_, HBHERecHitsHandle);
0058   if (!HBHERecHitsHandle.isValid()) {
0059     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
0060   } else {
0061     HBHERecHits = HBHERecHitsHandle.product();  // get a ptr to the product
0062   }
0063 
0064   iEvent.getByToken(tok_ho_, HORecHitsHandle);
0065   if (!HORecHitsHandle.isValid()) {
0066     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
0067   } else {
0068     HORecHits = HORecHitsHandle.product();  // get a ptr to the product
0069   }
0070 
0071   iEvent.getByToken(tok_hf_, HFRecHitsHandle);
0072   if (!HFRecHitsHandle.isValid()) {
0073     LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
0074   } else {
0075     HFRecHits = HFRecHitsHandle.product();  // get a ptr to the product
0076   }
0077 
0078   //Create empty output collections
0079   auto RecalibHBHERecHitCollection = std::make_unique<HBHERecHitCollection>();
0080   auto RecalibHFRecHitCollection = std::make_unique<HFRecHitCollection>();
0081   auto RecalibHORecHitCollection = std::make_unique<HORecHitCollection>();
0082 
0083   if (HBHERecHits) {
0084     HBHERecHitCollection::const_iterator itHBHE;
0085     for (itHBHE = HBHERecHits->begin(); itHBHE != HBHERecHits->end(); ++itHBHE) {
0086       // make the rechit with rescaled energy and put in the output collection
0087       float icalconst = (mapHcal_.get().find(itHBHE->id().rawId()))->second;
0088       icalconst = refactor_mean_ +
0089                   (icalconst - refactor_mean_) * refactor_;  //apply additional scaling factor (works if gaussian)
0090       HBHERecHit aHit(itHBHE->id(), itHBHE->energy() * icalconst, itHBHE->time());
0091 
0092       RecalibHBHERecHitCollection->push_back(aHit);
0093     }
0094   }
0095 
0096   if (HFRecHits) {
0097     HFRecHitCollection::const_iterator itHF;
0098     for (itHF = HFRecHits->begin(); itHF != HFRecHits->end(); ++itHF) {
0099       // make the rechit with rescaled energy and put in the output collection
0100       float icalconst = (mapHcal_.get().find(itHF->id().rawId()))->second;
0101       icalconst = refactor_mean_ +
0102                   (icalconst - refactor_mean_) * refactor_;  //apply additional scaling factor (works if gaussian)
0103       HFRecHit aHit(itHF->id(), itHF->energy() * icalconst, itHF->time());
0104 
0105       RecalibHFRecHitCollection->push_back(aHit);
0106     }
0107   }
0108 
0109   if (HORecHits) {
0110     HORecHitCollection::const_iterator itHO;
0111     for (itHO = HORecHits->begin(); itHO != HORecHits->end(); ++itHO) {
0112       // make the rechit with rescaled energy and put in the output collection
0113       float icalconst = (mapHcal_.get().find(itHO->id().rawId()))->second;
0114       icalconst = refactor_mean_ +
0115                   (icalconst - refactor_mean_) * refactor_;  //apply additional scaling factor (works if gaussian)
0116       HORecHit aHit(itHO->id(), itHO->energy() * icalconst, itHO->time());
0117 
0118       RecalibHORecHitCollection->push_back(aHit);
0119     }
0120   }
0121 
0122   //Put Recalibrated rechit in the event
0123   iEvent.put(std::move(RecalibHBHERecHitCollection), recalibHBHEHits_);
0124   iEvent.put(std::move(RecalibHFRecHitCollection), recalibHFHits_);
0125   iEvent.put(std::move(RecalibHORecHitCollection), recalibHOHits_);
0126 }