Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "CommonTools/CandAlgos/interface/ModifyObjectValueBase.h"
0002 
0003 #include "FWCore/Utilities/interface/InputTag.h"
0004 #include "FWCore/Utilities/interface/EDGetToken.h"
0005 
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/ESHandle.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 
0014 #include "DataFormats/Common/interface/ValueMap.h"
0015 
0016 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0017 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0018 
0019 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0020 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0021 
0022 #include "RecoEgamma/EgammaTools/interface/EgammaRegressionContainer.h"
0023 #include "RecoEgamma/EgammaTools/interface/EpCombinationTool.h"
0024 #include "RecoEgamma/EgammaTools/interface/EcalClusterLocal.h"
0025 
0026 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0027 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0028 
0029 /*
0030  * EGRegressionModifierDRN
0031  *
0032  * Object modifier to apply DRN regression.
0033  * Designed to be a drop-in replacement for EGRegressionModifierVX
0034  * 
0035  * Requires the appropriate DRNCorrectionProducerX(s) to also be in the path
0036  * You can specify which of reco::GsfElectron, reco::Photon, pat::Electron, pat::Photon 
0037  *    to apply corrections to in the config
0038  *
0039  */
0040 
0041 class EGRegressionModifierDRN : public ModifyObjectValueBase {
0042 public:
0043   EGRegressionModifierDRN(const edm::ParameterSet& conf, edm::ConsumesCollector& cc);
0044 
0045   void setEvent(const edm::Event&) final;
0046   void setEventContent(const edm::EventSetup&) final;
0047 
0048   void modifyObject(reco::GsfElectron&) const final;
0049   void modifyObject(reco::Photon&) const final;
0050 
0051   void modifyObject(pat::Electron&) const final;
0052   void modifyObject(pat::Photon&) const final;
0053 
0054 private:
0055   template <typename T>
0056   struct partVars {
0057     edm::InputTag source;
0058     edm::EDGetTokenT<edm::View<T>> token;
0059     const edm::View<T>* particles;
0060 
0061     edm::InputTag correctionsSource;
0062     edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> correctionsToken;
0063     const edm::ValueMap<std::pair<float, float>>* corrections;
0064 
0065     bool userFloat;
0066     std::string energyFloat, resFloat;
0067 
0068     unsigned i;
0069 
0070     partVars(const edm::ParameterSet& config, edm::ConsumesCollector& cc) {
0071       source = config.getParameter<edm::InputTag>("source");
0072       token = cc.consumes(source);
0073 
0074       correctionsSource = config.getParameter<edm::InputTag>("correctionsSource");
0075       correctionsToken = cc.consumes(correctionsSource);
0076 
0077       if (config.exists("energyFloat")) {
0078         userFloat = true;
0079         energyFloat = config.getParameter<std::string>("energyFloat");
0080         resFloat = config.getParameter<std::string>("resFloat");
0081       } else {
0082         userFloat = false;
0083       }
0084 
0085       i = 0;
0086     }
0087 
0088     const std::pair<float, float> getCorrection(T& part);
0089 
0090     const void doUserFloat(T& part, const std::pair<float, float>& correction) const {
0091       part.addUserFloat(energyFloat, correction.first);
0092       part.addUserFloat(resFloat, correction.second);
0093     }
0094   };
0095 
0096   std::unique_ptr<partVars<pat::Photon>> patPhotons_;
0097   std::unique_ptr<partVars<pat::Electron>> patElectrons_;
0098   std::unique_ptr<partVars<reco::Photon>> gedPhotons_;
0099   std::unique_ptr<partVars<reco::GsfElectron>> gsfElectrons_;
0100 };
0101 
0102 EGRegressionModifierDRN::EGRegressionModifierDRN(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0103     : ModifyObjectValueBase(conf) {
0104   if (conf.exists("patPhotons")) {
0105     patPhotons_ = std::make_unique<partVars<pat::Photon>>(conf.getParameterSet("patPhotons"), cc);
0106   }
0107 
0108   if (conf.exists("gedPhotons")) {
0109     gedPhotons_ = std::make_unique<partVars<reco::Photon>>(conf.getParameterSet("gedPhotons"), cc);
0110   }
0111 
0112   if (conf.exists("patElectrons")) {
0113     patElectrons_ = std::make_unique<partVars<pat::Electron>>(conf.getParameterSet("patElectrons"), cc);
0114   }
0115 
0116   if (conf.exists("gsfElectrons")) {
0117     gsfElectrons_ = std::make_unique<partVars<reco::GsfElectron>>(conf.getParameterSet("gsfElectrons"), cc);
0118   }
0119 }
0120 
0121 void EGRegressionModifierDRN::setEvent(const edm::Event& evt) {
0122   if (patElectrons_) {
0123     patElectrons_->particles = &evt.get(patElectrons_->token);
0124     patElectrons_->corrections = &evt.get(patElectrons_->correctionsToken);
0125     patElectrons_->i = 0;
0126   }
0127 
0128   if (patPhotons_) {
0129     patPhotons_->particles = &evt.get(patPhotons_->token);
0130     patPhotons_->corrections = &evt.get(patPhotons_->correctionsToken);
0131     patPhotons_->i = 0;
0132   }
0133 
0134   if (gsfElectrons_) {
0135     gsfElectrons_->particles = &evt.get(gsfElectrons_->token);
0136     gsfElectrons_->corrections = &evt.get(gsfElectrons_->correctionsToken);
0137     gsfElectrons_->i = 0;
0138   }
0139 
0140   if (gedPhotons_) {
0141     gedPhotons_->particles = &evt.get(gedPhotons_->token);
0142     gedPhotons_->corrections = &evt.get(gedPhotons_->correctionsToken);
0143     gedPhotons_->i = 0;
0144   }
0145 }
0146 
0147 void EGRegressionModifierDRN::setEventContent(const edm::EventSetup& iSetup) {}
0148 
0149 void EGRegressionModifierDRN::modifyObject(reco::GsfElectron& ele) const {
0150   if (!gsfElectrons_)
0151     return;
0152 
0153   const std::pair<float, float>& correction = gsfElectrons_->getCorrection(ele);
0154 
0155   if (correction.first > 0 && correction.second > 0) {
0156     ele.setCorrectedEcalEnergy(correction.first, true);
0157     ele.setCorrectedEcalEnergyError(correction.second);
0158   }
0159 
0160   throw cms::Exception("EGRegressionModifierDRN")
0161       << "Electron energy corrections not fully implemented yet:" << std::endl
0162       << "Still need E/p combination" << std::endl
0163       << "Do not enable DRN for electrons" << std::endl;
0164 }
0165 
0166 void EGRegressionModifierDRN::modifyObject(pat::Electron& ele) const {
0167   if (!patElectrons_)
0168     return;
0169 
0170   const std::pair<float, float>& correction = patElectrons_->getCorrection(ele);
0171 
0172   if (patElectrons_->userFloat) {
0173     patElectrons_->doUserFloat(ele, correction);
0174   } else if (correction.first > 0 && correction.second > 0) {
0175     ele.setCorrectedEcalEnergy(correction.first, true);
0176     ele.setCorrectedEcalEnergyError(correction.second);
0177   }
0178 
0179   throw cms::Exception("EGRegressionModifierDRN")
0180       << "Electron energy corrections not fully implemented yet:" << std::endl
0181       << "Still need E/p combination" << std::endl
0182       << "Do not enable DRN for electrons" << std::endl;
0183 }
0184 
0185 void EGRegressionModifierDRN::modifyObject(pat::Photon& pho) const {
0186   if (!patPhotons_)
0187     return;
0188 
0189   const std::pair<float, float>& correction = patPhotons_->getCorrection(pho);
0190 
0191   if (patPhotons_->userFloat) {
0192     patPhotons_->doUserFloat(pho, correction);
0193   } else if (correction.first > 0 && correction.second > 0) {
0194     pho.setCorrectedEnergy(pat::Photon::P4type::regression2, correction.first, correction.second, true);
0195   }
0196 }
0197 
0198 void EGRegressionModifierDRN::modifyObject(reco::Photon& pho) const {
0199   if (!gedPhotons_)
0200     return;
0201 
0202   const std::pair<float, float>& correction = gedPhotons_->getCorrection(pho);
0203 
0204   if (correction.first > 0 && correction.second > 0) {
0205     pho.setCorrectedEnergy(reco::Photon::P4type::regression2, correction.first, correction.second, true);
0206   }
0207 };
0208 
0209 template <typename T>
0210 const std::pair<float, float> EGRegressionModifierDRN::partVars<T>::getCorrection(T& part) {
0211   edm::Ptr<T> ptr = particles->ptrAt(i++);
0212 
0213   std::pair<float, float> correction = (*corrections)[ptr];
0214 
0215   return correction;
0216 }
0217 
0218 DEFINE_EDM_PLUGIN(ModifyObjectValueFactory, EGRegressionModifierDRN, "EGRegressionModifierDRN");