File indexing completed on 2024-04-06 12:19:24
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <vector>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/stream/EDProducer.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0032 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0033
0034 #include "FWCore/Utilities/interface/EDGetToken.h"
0035 #include "FWCore/Utilities/interface/InputTag.h"
0036
0037 #include "DataFormats/Common/interface/Handle.h"
0038
0039 #include "JetMETCorrections/JetCorrector/interface/JetCorrector.h"
0040 #include "JetMETCorrections/JetCorrector/interface/JetCorrectorImpl.h"
0041
0042
0043
0044
0045 namespace {
0046 class ChainedJetCorrectorImpl : public reco::JetCorrectorImpl {
0047 public:
0048 ChainedJetCorrectorImpl(std::vector<reco::JetCorrector const*> correctors)
0049 : jetCorrectors_(std::move(correctors)) {}
0050
0051 double correction(const LorentzVector& fJet) const override;
0052
0053
0054 double correction(const reco::Jet& fJet) const override;
0055
0056
0057 double correction(const reco::Jet& fJet, const edm::RefToBase<reco::Jet>& fJetRef) const override;
0058
0059
0060 bool refRequired() const override;
0061
0062 private:
0063 std::vector<reco::JetCorrector const*> jetCorrectors_;
0064 };
0065
0066 double ChainedJetCorrectorImpl::correction(const LorentzVector& fJet) const {
0067 LorentzVector jet = fJet;
0068 double result = 1;
0069 for (auto cor : jetCorrectors_) {
0070 double scale = cor->correction(jet);
0071 jet *= scale;
0072 result *= scale;
0073 }
0074 return result;
0075 }
0076
0077
0078 double ChainedJetCorrectorImpl::correction(const reco::Jet& fJet) const {
0079 std::unique_ptr<reco::Jet> jet(dynamic_cast<reco::Jet*>(fJet.clone()));
0080 double result = 1;
0081 for (auto cor : jetCorrectors_) {
0082 double scale = cor->correction(*jet);
0083 jet->scaleEnergy(scale);
0084 result *= scale;
0085 }
0086 return result;
0087 }
0088
0089
0090 double ChainedJetCorrectorImpl::correction(const reco::Jet& fJet, const edm::RefToBase<reco::Jet>& fJetRef) const {
0091 std::unique_ptr<reco::Jet> jet(dynamic_cast<reco::Jet*>(fJet.clone()));
0092 double result = 1;
0093 for (auto cor : jetCorrectors_) {
0094 double scale = cor->correction(*jet, fJetRef);
0095 jet->scaleEnergy(scale);
0096 result *= scale;
0097 }
0098 return result;
0099 }
0100
0101
0102 bool ChainedJetCorrectorImpl::refRequired() const {
0103 for (auto cor : jetCorrectors_) {
0104 if (cor->refRequired())
0105 return true;
0106 }
0107 return false;
0108 }
0109
0110 }
0111
0112 class ChainedJetCorrectorProducer : public edm::stream::EDProducer<> {
0113 public:
0114 explicit ChainedJetCorrectorProducer(const edm::ParameterSet&);
0115
0116 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0117
0118 private:
0119 void produce(edm::Event&, const edm::EventSetup&) override;
0120
0121
0122 std::vector<edm::EDGetTokenT<reco::JetCorrector>> correctorTokens_;
0123 };
0124
0125
0126
0127
0128
0129
0130
0131
0132
0133
0134
0135
0136 ChainedJetCorrectorProducer::ChainedJetCorrectorProducer(const edm::ParameterSet& iConfig) {
0137
0138 produces<reco::JetCorrector>();
0139
0140 auto const& tags = iConfig.getParameter<std::vector<edm::InputTag>>("correctors");
0141 correctorTokens_.reserve(tags.size());
0142
0143 for (auto const& tag : tags) {
0144 correctorTokens_.emplace_back(consumes<reco::JetCorrector>(tag));
0145 }
0146 }
0147
0148
0149
0150
0151
0152
0153 void ChainedJetCorrectorProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0154 using namespace edm;
0155
0156 std::vector<reco::JetCorrector const*> correctors;
0157 correctors.reserve(correctorTokens_.size());
0158
0159 for (auto const& token : correctorTokens_) {
0160 edm::Handle<reco::JetCorrector> hCorrector;
0161 iEvent.getByToken(token, hCorrector);
0162 correctors.emplace_back(&(*hCorrector));
0163 }
0164
0165 std::unique_ptr<reco::JetCorrector> pCorr(new reco::JetCorrector(
0166 std::unique_ptr<reco::JetCorrectorImpl>(new ChainedJetCorrectorImpl(std::move(correctors)))));
0167 iEvent.put(std::move(pCorr));
0168 }
0169
0170
0171 void ChainedJetCorrectorProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0172
0173
0174 edm::ParameterSetDescription desc;
0175 desc.add<std::vector<edm::InputTag>>("correctors");
0176 descriptions.addDefault(desc);
0177 }
0178
0179
0180 DEFINE_FWK_MODULE(ChainedJetCorrectorProducer);