Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-12 23:42:04

0001 /** \class EcalRecalibRecHitProducer
0002  *   produce ECAL rechits from uncalibrated rechits
0003  *
0004  *  \author Federico Ferri, University of Milano Bicocca and INFN
0005  *
0006  **/
0007 
0008 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbRecord.h"
0009 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbService.h"
0010 #include "CondFormats/DataRecord/interface/EcalADCToGeVConstantRcd.h"
0011 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
0012 #include "CondFormats/EcalObjects/interface/EcalADCToGeVConstant.h"
0013 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0016 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0017 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
0018 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0019 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
0020 #include "FWCore/Framework/interface/ESHandle.h"
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/global/EDProducer.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0026 #include "FWCore/Utilities/interface/ESGetToken.h"
0027 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalRecHitAbsAlgo.h"
0028 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalRecHitSimpleAlgo.h"
0029 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0030 
0031 #include <cmath>
0032 #include <iostream>
0033 #include <vector>
0034 
0035 class EcalRecalibRecHitProducer : public edm::global::EDProducer<> {
0036 public:
0037   explicit EcalRecalibRecHitProducer(const edm::ParameterSet& ps);
0038   void produce(edm::StreamID sid, edm::Event& evt, const edm::EventSetup& es) const override;
0039 
0040   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0041 
0042 private:
0043   const edm::InputTag EBRecHitCollection_;
0044   const edm::InputTag EERecHitCollection_;
0045   const edm::EDGetTokenT<EBRecHitCollection> EBRecHitToken_;
0046   const edm::EDGetTokenT<EERecHitCollection> EERecHitToken_;
0047 
0048   const std::string EBRecalibRecHitCollection_;  // secondary name to be given to EB collection of hits
0049   const std::string EERecalibRecHitCollection_;  // secondary name to be given to EE collection of hits
0050 
0051   const bool doEnergyScale_;
0052   const bool doIntercalib_;
0053   const bool doLaserCorrections_;
0054   const bool doEnergyScaleInverse_;
0055   const bool doIntercalibInverse_;
0056   const bool doLaserCorrectionsInverse_;
0057 
0058   edm::ESGetToken<EcalADCToGeVConstant, EcalADCToGeVConstantRcd> ecalADCToGeVConstantToken_;
0059   edm::ESGetToken<EcalIntercalibConstants, EcalIntercalibConstantsRcd> ecalIntercalibConstantsToken_;
0060   edm::ESGetToken<EcalLaserDbService, EcalLaserDbRecord> ecalLaserDBServiceToken_;
0061 };
0062 
0063 EcalRecalibRecHitProducer::EcalRecalibRecHitProducer(const edm::ParameterSet& ps)
0064     : EBRecHitCollection_(ps.getParameter<edm::InputTag>("EBRecHitCollection")),
0065       EERecHitCollection_(ps.getParameter<edm::InputTag>("EERecHitCollection")),
0066       EBRecHitToken_((not EBRecHitCollection_.label().empty()) ? consumes<EBRecHitCollection>(EBRecHitCollection_)
0067                                                                : edm::EDGetTokenT<EBRecHitCollection>()),
0068       EERecHitToken_((not EERecHitCollection_.label().empty()) ? consumes<EERecHitCollection>(EERecHitCollection_)
0069                                                                : edm::EDGetTokenT<EERecHitCollection>()),
0070       EBRecalibRecHitCollection_(ps.getParameter<std::string>("EBRecalibRecHitCollection")),
0071       EERecalibRecHitCollection_(ps.getParameter<std::string>("EERecalibRecHitCollection")),
0072       doEnergyScale_(ps.getParameter<bool>("doEnergyScale")),
0073       doIntercalib_(ps.getParameter<bool>("doIntercalib")),
0074       doLaserCorrections_(ps.getParameter<bool>("doLaserCorrections")),
0075 
0076       doEnergyScaleInverse_(ps.getParameter<bool>("doEnergyScaleInverse")),
0077       doIntercalibInverse_(ps.getParameter<bool>("doIntercalibInverse")),
0078       doLaserCorrectionsInverse_(ps.getParameter<bool>("doLaserCorrectionsInverse")),
0079       ecalLaserDBServiceToken_(esConsumes<EcalLaserDbService, EcalLaserDbRecord>()) {
0080   if (doEnergyScale_) {
0081     ecalADCToGeVConstantToken_ = esConsumes<EcalADCToGeVConstant, EcalADCToGeVConstantRcd>();
0082   }
0083   if (doIntercalib_) {
0084     ecalIntercalibConstantsToken_ = esConsumes<EcalIntercalibConstants, EcalIntercalibConstantsRcd>();
0085   }
0086   produces<EBRecHitCollection>(EBRecalibRecHitCollection_);
0087   produces<EERecHitCollection>(EERecalibRecHitCollection_);
0088 }
0089 
0090 void EcalRecalibRecHitProducer::produce(edm::StreamID sid, edm::Event& evt, const edm::EventSetup& es) const {
0091   using namespace edm;
0092   Handle<EBRecHitCollection> pEBRecHits;
0093   Handle<EERecHitCollection> pEERecHits;
0094 
0095   const EBRecHitCollection* EBRecHits = nullptr;
0096   const EERecHitCollection* EERecHits = nullptr;
0097 
0098   if (not EBRecHitCollection_.label().empty()) {
0099     evt.getByToken(EBRecHitToken_, pEBRecHits);
0100     EBRecHits = pEBRecHits.product();  // get a ptr to the product
0101   }
0102   if (not EERecHitCollection_.label().empty()) {
0103     evt.getByToken(EERecHitToken_, pEERecHits);
0104     EERecHits = pEERecHits.product();  // get a ptr to the product
0105   }
0106 
0107   // collection of rechits to put in the event
0108   auto EBRecalibRecHits = std::make_unique<EBRecHitCollection>();
0109   auto EERecalibRecHits = std::make_unique<EERecHitCollection>();
0110 
0111   // now fetch all conditions we need to make rechits
0112   // ADC to GeV constant
0113   edm::ESHandle<EcalADCToGeVConstant> pAgc;
0114   const EcalADCToGeVConstant* agc = nullptr;
0115   float agc_eb = 1.;
0116   float agc_ee = 1.;
0117   if (doEnergyScale_) {
0118     pAgc = es.getHandle(ecalADCToGeVConstantToken_);
0119     agc = pAgc.product();
0120     // use this value in the algorithm
0121     agc_eb = float(agc->getEBValue());
0122     agc_ee = float(agc->getEEValue());
0123   }
0124   // Intercalib constants
0125   edm::ESHandle<EcalIntercalibConstants> pIcal;
0126   const EcalIntercalibConstants* ical = nullptr;
0127   if (doIntercalib_) {
0128     pIcal = es.getHandle(ecalIntercalibConstantsToken_);
0129     ical = pIcal.product();
0130   }
0131   // Laser corrections
0132   edm::ESHandle<EcalLaserDbService> pLaser = es.getHandle(ecalLaserDBServiceToken_);
0133 
0134   if (doEnergyScaleInverse_) {
0135     agc_eb = 1.0 / agc_eb;
0136     agc_ee = 1.0 / agc_ee;
0137   }
0138 
0139   if (EBRecHits) {
0140     // loop over uncalibrated rechits to make calibrated ones
0141     for (EBRecHitCollection::const_iterator it = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
0142       EcalIntercalibConstant icalconst = 1.;
0143       if (doIntercalib_) {
0144         // find intercalib constant for this xtal
0145         const EcalIntercalibConstantMap& icalMap = ical->getMap();
0146         EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
0147         if (icalit != icalMap.end()) {
0148           icalconst = (*icalit);
0149         } else {
0150           edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(it->id())
0151                                            << "! something wrong with EcalIntercalibConstants in your DB? ";
0152         }
0153       }
0154       // get laser coefficient
0155       float lasercalib = 1;
0156       if (doLaserCorrections_) {
0157         lasercalib = pLaser->getLaserCorrection(EBDetId(it->id()), evt.time());
0158       }
0159 
0160       // make the rechit and put in the output collection
0161       // must implement op= for EcalRecHit
0162 
0163       if (doIntercalibInverse_) {
0164         icalconst = 1.0 / icalconst;
0165       }
0166       if (doLaserCorrectionsInverse_) {
0167         lasercalib = 1.0 / lasercalib;
0168       }
0169 
0170       EcalRecHit aHit((*it).id(), (*it).energy() * agc_eb * icalconst * lasercalib, (*it).time());
0171       EBRecalibRecHits->push_back(aHit);
0172     }
0173   }
0174 
0175   if (EERecHits) {
0176     // loop over uncalibrated rechits to make calibrated ones
0177     for (EERecHitCollection::const_iterator it = EERecHits->begin(); it != EERecHits->end(); ++it) {
0178       // find intercalib constant for this xtal
0179       EcalIntercalibConstant icalconst = 1.;
0180       if (doIntercalib_) {
0181         const EcalIntercalibConstantMap& icalMap = ical->getMap();
0182         EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
0183         if (icalit != icalMap.end()) {
0184           icalconst = (*icalit);
0185         } else {
0186           edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EEDetId(it->id())
0187                                            << "! something wrong with EcalIntercalibConstants in your DB? ";
0188         }
0189       }
0190       // get laser coefficient
0191       float lasercalib = 1;
0192       if (doLaserCorrections_) {
0193         lasercalib = pLaser->getLaserCorrection(EEDetId(it->id()), evt.time());
0194       }
0195 
0196       if (doIntercalibInverse_) {
0197         icalconst = 1.0 / icalconst;
0198       }
0199       if (doLaserCorrectionsInverse_) {
0200         lasercalib = 1.0 / lasercalib;
0201       }
0202 
0203       // make the rechit and put in the output collection
0204       EcalRecHit aHit((*it).id(), (*it).energy() * agc_ee * icalconst * lasercalib, (*it).time());
0205       EERecalibRecHits->push_back(aHit);
0206     }
0207   }
0208   // put the collection of recunstructed hits in the event
0209   LogInfo("EcalRecalibRecHitInfo") << "total # EB re-calibrated rechits: " << EBRecalibRecHits->size();
0210   LogInfo("EcalRecalibRecHitInfo") << "total # EE re-calibrated rechits: " << EERecalibRecHits->size();
0211 
0212   evt.put(std::move(EBRecalibRecHits), EBRecalibRecHitCollection_);
0213   evt.put(std::move(EERecalibRecHits), EERecalibRecHitCollection_);
0214 }
0215 
0216 void EcalRecalibRecHitProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0217   edm::ParameterSetDescription desc;
0218   desc.add<edm::InputTag>("EBRecHitCollection", edm::InputTag("ecalRecHit", "EcalRecHitsEB"));
0219   desc.add<edm::InputTag>("EERecHitCollection", edm::InputTag("ecalRecHit", "EcalRecHitsEE"));
0220   desc.add<std::string>("EBRecalibRecHitCollection", "EcalRecHitsEB");
0221   desc.add<std::string>("EERecalibRecHitCollection", "EcalRecHitsEE");
0222   desc.add<bool>("doEnergyScale", false);
0223   desc.add<bool>("doIntercalib", false);
0224   desc.add<bool>("doLaserCorrections", false);
0225   desc.add<bool>("doEnergyScaleInverse", false);
0226   desc.add<bool>("doIntercalibInverse", false);
0227   desc.add<bool>("doLaserCorrectionsInverse", false);
0228   descriptions.addWithDefaultLabel(desc);
0229 }
0230 
0231 #include "FWCore/Framework/interface/MakerMacros.h"
0232 DEFINE_FWK_MODULE(EcalRecalibRecHitProducer);