Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:55

0001 /**
0002   \class    pat::TauJetCorrFactorsProducer TauJetCorrFactorsProducer.h "PhysicsTools/PatAlgos/interface/TauJetCorrFactorsProducer.h"
0003   \brief    Produces a ValueMap between TauJetCorrFactors and the originating reco taus
0004 
0005    The TauJetCorrFactorsProducer produces a set of tau-jet energy correction factors, defined in the class pat::TauJetCorrFactors.
0006    This vector is linked to the originating reco taus through an edm::ValueMap. The initializing parameters of the module can be found
0007    in the recoLayer1/tauJetCorrFactors_cfi.py of the PatAlgos package. In the standard PAT workflow the module has to be run
0008    before the creation of the pat::Tau. The edm::ValueMap will then be embedded into the pat::Tau.
0009 
0010    Jets corrected up to a given correction level can then be accessed via the pat::Tau member function correctedJet. For
0011    more details have a look into the class description of the pat::Tau.
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     /// value map for JetCorrFactors (to be written into the event)
0040     typedef edm::ValueMap<pat::TauJetCorrFactors> JetCorrFactorsMap;
0041 
0042   public:
0043     /// default constructor
0044     explicit TauJetCorrFactorsProducer(const edm::ParameterSet&);
0045     /// default destructor
0046     ~TauJetCorrFactorsProducer() override{};
0047 
0048     /// everything that needs to be done per event
0049     void produce(edm::Event&, const edm::EventSetup&) override;
0050 
0051   private:
0052     /// return the jec parameters as input to the FactorizedJetCorrector for different flavors
0053     std::vector<JetCorrectorParameters> params(const JetCorrectorParametersCollection&,
0054                                                const std::vector<std::string>&) const;
0055 
0056     /// evaluate jet correction factor up to a given level
0057     float evaluate(edm::View<reco::BaseTau>::const_iterator&, std::shared_ptr<FactorizedJetCorrector>&, int);
0058 
0059   private:
0060     /// python label of this TauJetCorrFactorsProducer module
0061     std::string moduleLabel_;
0062 
0063     /// input tau-jet collection
0064     edm::EDGetTokenT<edm::View<reco::BaseTau> > srcToken_;
0065 
0066     /// mapping of reconstructed tau decay modes to payloads
0067     typedef std::vector<int> vint;
0068     struct payloadMappingType {
0069       /// reconstructed tau decay modes associated to this payload,
0070       /// as defined in DataFormats/TauReco/interface/PFTau.h
0071       vint decayModes_;
0072 
0073       /// payload label
0074       std::string payload_;
0075     };
0076     std::vector<payloadMappingType> payloadMappings_;
0077 
0078     /// payload to be used for decay modes not explicitely specified
0079     ///
0080     /// NOTE: no decay mode reconstruction implemented for CaloTaus so far
0081     ///      --> this payload is used for all CaloTaus
0082     ///
0083     std::string defaultPayload_;
0084 
0085     /// jec levels
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 }  // namespace pat
0093 
0094 /// value map for JetCorrFactors (to be written into the event)
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   // get tau-jet collection from the event
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   // fill the tauJetCorrFactors
0160   std::vector<TauJetCorrFactors> tauJetCorrections;
0161   for (edm::View<reco::BaseTau>::const_iterator tauJet = tauJets->begin(); tauJet != tauJets->end(); ++tauJet) {
0162     // the TauJetCorrFactors::CorrectionFactor is a std::pair<std::string, float>
0163     // the string corresponds to the label of the correction level, the float to the tau-jet energy correction factor.
0164     // The first correction level is predefined with label 'Uncorrected'. The correction factor is 1.
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     // retrieve JEC parameters from the DB and build a new corrector,
0190     // in case it does not exist already for current payload
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     // evaluate tau-jet energy corrections
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       // push back the set of JetCorrFactors: the first entry corresponds to the label
0206       // of the correction level. The second parameter corresponds to the jec factor.
0207       // In the default configuration the CorrectionFactor will look like this:
0208       //   'Uncorrected' : 1 ;
0209       //   'L2Relative'  : x ;
0210       //   'L3Absolute'  : x ;
0211       jec.push_back(std::make_pair(corrLevel.substr(0, corrLevel.find('_')), jecFactor));
0212     }
0213 
0214     // create the actual object with the scale factors we want the valuemap to refer to
0215     // moduleLabel_ corresponds to the python label of the TauJetCorrFactorsProducer module instance
0216     TauJetCorrFactors tauJetCorrection(moduleLabel_, jec);
0217     tauJetCorrections.push_back(tauJetCorrection);
0218   }
0219 
0220   // build the valuemap
0221   auto jecMapping = std::make_unique<TauJetCorrFactorsMap>();
0222   TauJetCorrFactorsMap::Filler filler(*jecMapping);
0223   // tauJets and tauJetCorrections vectors have their indices aligned by construction
0224   filler.insert(tauJets, tauJetCorrections.begin(), tauJetCorrections.end());
0225   filler.fill();  // do the actual filling
0226 
0227   // add valuemap to the event
0228   evt.put(std::move(jecMapping));
0229 }
0230 
0231 #include "FWCore/Framework/interface/MakerMacros.h"
0232 
0233 DEFINE_FWK_MODULE(TauJetCorrFactorsProducer);