Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-10-19 10:28:34

0001 #ifndef JetCorrectionProducer_h
0002 #define JetCorrectionProducer_h
0003 
0004 #include <sstream>
0005 #include <string>
0006 #include <vector>
0007 
0008 #include "CommonTools/Utils/interface/PtComparator.h"
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "DataFormats/Common/interface/Ref.h"
0011 #include "DataFormats/Common/interface/RefToBase.h"
0012 #include "FWCore/Framework/interface/stream/EDProducer.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "JetMETCorrections/Objects/interface/JetCorrectionsRecord.h"
0019 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
0020 
0021 namespace edm {
0022   class ParameterSet;
0023 }
0024 
0025 class JetCorrector;
0026 
0027 namespace cms {
0028   template <class T>
0029   class JetCorrectionProducer : public edm::stream::EDProducer<> {
0030   public:
0031     typedef std::vector<T> JetCollection;
0032     explicit JetCorrectionProducer(const edm::ParameterSet& fParameters);
0033     ~JetCorrectionProducer() override {}
0034     void produce(edm::Event&, const edm::EventSetup&) override;
0035 
0036   private:
0037     edm::EDGetTokenT<JetCollection> mInput;
0038     std::vector<std::string> mCorrectorNames;
0039     std::vector<edm::ESGetToken<JetCorrector, JetCorrectionsRecord>> mCorrectorTokens;
0040     // cache
0041     std::vector<const JetCorrector*> mCorrectors;
0042     unsigned long long mCacheId;
0043     bool mVerbose;
0044   };
0045 }  // namespace cms
0046 
0047 // ---------- implementation ----------
0048 
0049 namespace cms {
0050 
0051   template <class T>
0052   JetCorrectionProducer<T>::JetCorrectionProducer(const edm::ParameterSet& fConfig)
0053       : mInput(consumes<JetCollection>(fConfig.getParameter<edm::InputTag>("src"))),
0054         mCorrectorNames(fConfig.getParameter<std::vector<std::string>>("correctors")),
0055         mCorrectors(mCorrectorNames.size(), nullptr),
0056         mCacheId(0),
0057         mVerbose(fConfig.getUntrackedParameter<bool>("verbose", false)) {
0058     mCorrectorTokens.reserve(mCorrectorNames.size());
0059     for (auto const& n : mCorrectorNames) {
0060       mCorrectorTokens.emplace_back(esConsumes(edm::ESInputTag("", n)));
0061     }
0062     std::string alias = fConfig.getUntrackedParameter<std::string>("alias", "");
0063     if (alias.empty())
0064       produces<JetCollection>();
0065     else
0066       produces<JetCollection>().setBranchAlias(alias);
0067   }
0068 
0069   template <class T>
0070   void JetCorrectionProducer<T>::produce(edm::Event& fEvent, const edm::EventSetup& fSetup) {
0071     // look for correctors
0072     const JetCorrectionsRecord& record = fSetup.get<JetCorrectionsRecord>();
0073     if (record.cacheIdentifier() != mCacheId) {  // need to renew cache
0074       for (unsigned i = 0; i < mCorrectorTokens.size(); i++) {
0075         mCorrectors[i] = &record.get(mCorrectorTokens[i]);
0076       }
0077       mCacheId = record.cacheIdentifier();
0078     }
0079     edm::Handle<JetCollection> jets;                           //Define Inputs
0080     fEvent.getByToken(mInput, jets);                           //Get Inputs
0081     std::unique_ptr<JetCollection> result(new JetCollection);  //Corrected jets
0082     typename JetCollection::const_iterator jet;
0083     for (jet = jets->begin(); jet != jets->end(); jet++) {
0084       const T* referenceJet = &*jet;
0085       int index = jet - jets->begin();
0086       edm::RefToBase<reco::Jet> jetRef(edm::Ref<JetCollection>(jets, index));
0087       T correctedJet = *jet;  //copy original jet
0088       if (mVerbose)
0089         std::cout << "JetCorrectionProducer::produce-> original jet: " << jet->print() << std::endl;
0090       for (unsigned i = 0; i < mCorrectors.size(); ++i) {
0091         if (!(mCorrectors[i]->vectorialCorrection())) {
0092           // Scalar correction
0093           double scale = 1.;
0094           if (!(mCorrectors[i]->refRequired()))
0095             scale = mCorrectors[i]->correction(*referenceJet, fEvent, fSetup);
0096           else
0097             scale = mCorrectors[i]->correction(*referenceJet, jetRef, fEvent, fSetup);
0098           if (mVerbose)
0099             std::cout << "JetCorrectionProducer::produce-> Corrector # " << i << ", correction factor: " << scale
0100                       << std::endl;
0101           correctedJet.scaleEnergy(scale);  // apply scalar correction
0102           referenceJet = &correctedJet;
0103         } else {
0104           // Vectorial correction
0105           JetCorrector::LorentzVector corrected;
0106           double scale = mCorrectors[i]->correction(*referenceJet, jetRef, fEvent, fSetup, corrected);
0107           if (mVerbose)
0108             std::cout << "JetCorrectionProducer::produce-> Corrector # " << i << ", correction factor: " << scale
0109                       << std::endl;
0110           correctedJet.setP4(corrected);  // apply vectorial correction
0111           referenceJet = &correctedJet;
0112         }
0113       }
0114       if (mVerbose)
0115         std::cout << "JetCorrectionProducer::produce-> corrected jet: " << correctedJet.print() << std::endl;
0116       result->push_back(correctedJet);
0117     }
0118     NumericSafeGreaterByPt<T> compJets;
0119     // reorder corrected jets
0120     std::sort(result->begin(), result->end(), compJets);
0121     // put corrected jet collection into event
0122     fEvent.put(std::move(result));
0123   }
0124 
0125 }  // namespace cms
0126 
0127 #endif