File indexing completed on 2024-09-07 04:37:22
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "JetMETCorrections/Objects/interface/JetCorrectionsRecord.h"
0016 #include "FWCore/Framework/interface/ESHandle.h"
0017 #include "DataFormats/TauReco/interface/PFTau.h"
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/Framework/interface/EventSetup.h"
0020 #include "FWCore/Utilities/interface/InputTag.h"
0021 #include "FWCore/Framework/interface/stream/EDProducer.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "CondFormats/JetMETObjects/interface/FactorizedJetCorrector.h"
0024 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
0025 #include "DataFormats/Common/interface/View.h"
0026 #include "DataFormats/TauReco/interface/BaseTau.h"
0027 #include "DataFormats/PatCandidates/interface/TauJetCorrFactors.h"
0028 #include "DataFormats/Common/interface/ValueMap.h"
0029
0030 #include <map>
0031 #include <memory>
0032 #include <string>
0033 #include <unordered_map>
0034
0035 namespace pat {
0036
0037 class TauJetCorrFactorsProducer : public edm::stream::EDProducer<> {
0038 public:
0039
0040 typedef edm::ValueMap<pat::TauJetCorrFactors> JetCorrFactorsMap;
0041
0042 public:
0043
0044 explicit TauJetCorrFactorsProducer(const edm::ParameterSet&);
0045
0046 ~TauJetCorrFactorsProducer() override {}
0047
0048
0049 void produce(edm::Event&, const edm::EventSetup&) override;
0050
0051 private:
0052
0053 std::vector<JetCorrectorParameters> params(const JetCorrectorParametersCollection&,
0054 const std::vector<std::string>&) const;
0055
0056
0057 float evaluate(edm::View<reco::BaseTau>::const_iterator&, std::shared_ptr<FactorizedJetCorrector>&, int);
0058
0059 private:
0060
0061 std::string moduleLabel_;
0062
0063
0064 edm::EDGetTokenT<edm::View<reco::BaseTau> > srcToken_;
0065
0066
0067 typedef std::vector<int> vint;
0068 struct payloadMappingType {
0069
0070
0071 vint decayModes_;
0072
0073
0074 std::string payload_;
0075 };
0076 std::vector<payloadMappingType> payloadMappings_;
0077
0078
0079
0080
0081
0082
0083 std::string defaultPayload_;
0084
0085
0086 typedef std::vector<std::string> vstring;
0087 vstring levels_;
0088
0089 using CorrectionToken = edm::ESGetToken<JetCorrectorParametersCollection, JetCorrectionsRecord>;
0090 std::unordered_map<std::string, CorrectionToken> payloadToTokens_;
0091 };
0092 }
0093
0094
0095 typedef edm::ValueMap<pat::TauJetCorrFactors> TauJetCorrFactorsMap;
0096
0097 using namespace pat;
0098
0099 TauJetCorrFactorsProducer::TauJetCorrFactorsProducer(const edm::ParameterSet& cfg)
0100 : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
0101 srcToken_(consumes<edm::View<reco::BaseTau> >(cfg.getParameter<edm::InputTag>("src"))),
0102 levels_(cfg.getParameter<std::vector<std::string> >("levels")) {
0103 typedef std::vector<edm::ParameterSet> vParameterSet;
0104 vParameterSet parameters = cfg.getParameter<vParameterSet>("parameters");
0105 for (vParameterSet::const_iterator param = parameters.begin(); param != parameters.end(); ++param) {
0106 payloadMappingType payloadMapping;
0107
0108 payloadMapping.payload_ = param->getParameter<std::string>("payload");
0109
0110 vstring decayModes_string = param->getParameter<vstring>("decayModes");
0111 for (vstring::const_iterator decayMode = decayModes_string.begin(); decayMode != decayModes_string.end();
0112 ++decayMode) {
0113 if ((*decayMode) == "*") {
0114 defaultPayload_ = payloadMapping.payload_;
0115 } else {
0116 payloadMapping.decayModes_.push_back(atoi(decayMode->data()));
0117 }
0118 }
0119
0120 if (!payloadMapping.decayModes_.empty()) {
0121 payloadMappings_.push_back(payloadMapping);
0122 payloadToTokens_.emplace(payloadMapping.payload_, CorrectionToken());
0123 }
0124 }
0125 for (auto& payloadToken : payloadToTokens_) {
0126 payloadToken.second = esConsumes(edm::ESInputTag("", payloadToken.first));
0127 }
0128
0129 produces<TauJetCorrFactorsMap>();
0130 }
0131
0132 std::vector<JetCorrectorParameters> TauJetCorrFactorsProducer::params(
0133 const JetCorrectorParametersCollection& jecParameters, const std::vector<std::string>& levels) const {
0134 std::vector<JetCorrectorParameters> retVal;
0135 for (std::vector<std::string>::const_iterator corrLevel = levels.begin(); corrLevel != levels.end(); ++corrLevel) {
0136 const JetCorrectorParameters& jecParameter_level = jecParameters[*corrLevel];
0137 retVal.push_back(jecParameter_level);
0138 }
0139 return retVal;
0140 }
0141
0142 float TauJetCorrFactorsProducer::evaluate(edm::View<reco::BaseTau>::const_iterator& tauJet,
0143 std::shared_ptr<FactorizedJetCorrector>& corrector,
0144 int corrLevel) {
0145 corrector->setJetEta(tauJet->eta());
0146 corrector->setJetPt(tauJet->pt());
0147 corrector->setJetE(tauJet->energy());
0148 return corrector->getSubCorrections()[corrLevel];
0149 }
0150
0151 void TauJetCorrFactorsProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0152
0153 edm::Handle<edm::View<reco::BaseTau> > tauJets;
0154 evt.getByToken(srcToken_, tauJets);
0155
0156 typedef std::shared_ptr<FactorizedJetCorrector> FactorizedJetCorrectorPtr;
0157 std::map<std::string, FactorizedJetCorrectorPtr> correctorMapping;
0158
0159
0160 std::vector<TauJetCorrFactors> tauJetCorrections;
0161 for (edm::View<reco::BaseTau>::const_iterator tauJet = tauJets->begin(); tauJet != tauJets->end(); ++tauJet) {
0162
0163
0164
0165 std::vector<TauJetCorrFactors::CorrectionFactor> jec;
0166 jec.push_back(std::make_pair(std::string("Uncorrected"), 1.0));
0167
0168 if (levels_.empty())
0169 throw cms::Exception("No JECFactors")
0170 << "You request to create a jetCorrFactors object with no JEC Levels indicated. \n"
0171 << "This makes no sense, either you should correct this or drop the module from \n"
0172 << "the sequence.";
0173
0174 std::string payload = defaultPayload_;
0175 if (dynamic_cast<const reco::PFTau*>(&(*tauJet))) {
0176 const reco::PFTau* pfTauJet = dynamic_cast<const reco::PFTau*>(&(*tauJet));
0177 for (std::vector<payloadMappingType>::const_iterator payloadMapping = payloadMappings_.begin();
0178 payloadMapping != payloadMappings_.end();
0179 ++payloadMapping) {
0180 for (vint::const_iterator decayMode = payloadMapping->decayModes_.begin();
0181 decayMode != payloadMapping->decayModes_.end();
0182 ++decayMode) {
0183 if (pfTauJet->decayMode() == (*decayMode))
0184 payload = payloadMapping->payload_;
0185 }
0186 }
0187 }
0188
0189
0190
0191 if (correctorMapping.find(payload) == correctorMapping.end()) {
0192 auto const& jecParameters = es.getData(payloadToTokens_[payload]);
0193
0194 correctorMapping[payload] = std::make_shared<FactorizedJetCorrector>(params(jecParameters, levels_));
0195 }
0196 FactorizedJetCorrectorPtr& corrector = correctorMapping[payload];
0197
0198
0199 size_t numLevels = levels_.size();
0200 for (size_t idx = 0; idx < numLevels; ++idx) {
0201 const std::string& corrLevel = levels_[idx];
0202
0203 float jecFactor = evaluate(tauJet, corrector, idx);
0204
0205
0206
0207
0208
0209
0210
0211 jec.push_back(std::make_pair(corrLevel.substr(0, corrLevel.find('_')), jecFactor));
0212 }
0213
0214
0215
0216 TauJetCorrFactors tauJetCorrection(moduleLabel_, jec);
0217 tauJetCorrections.push_back(tauJetCorrection);
0218 }
0219
0220
0221 auto jecMapping = std::make_unique<TauJetCorrFactorsMap>();
0222 TauJetCorrFactorsMap::Filler filler(*jecMapping);
0223
0224 filler.insert(tauJets, tauJetCorrections.begin(), tauJetCorrections.end());
0225 filler.fill();
0226
0227
0228 evt.put(std::move(jecMapping));
0229 }
0230
0231 #include "FWCore/Framework/interface/MakerMacros.h"
0232
0233 DEFINE_FWK_MODULE(TauJetCorrFactorsProducer);