Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:43

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