Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:26

0001 #include "JetMETCorrections/Type1MET/plugins/PFCandMETcorrInputProducer.h"
0002 
0003 #include <memory>
0004 
0005 #include "DataFormats/Common/interface/View.h"
0006 
0007 PFCandMETcorrInputProducer::PFCandMETcorrInputProducer(const edm::ParameterSet& cfg)
0008     : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
0009   token_ = consumes<edm::View<reco::Candidate>>(cfg.getParameter<edm::InputTag>("src"));
0010   edm::InputTag srcWeights = cfg.getParameter<edm::InputTag>("srcWeights");
0011   if (!srcWeights.label().empty())
0012     weightsToken_ = consumes<edm::ValueMap<float>>(srcWeights);
0013 
0014   if (cfg.exists("binning")) {
0015     typedef std::vector<edm::ParameterSet> vParameterSet;
0016     vParameterSet cfgBinning = cfg.getParameter<vParameterSet>("binning");
0017     for (vParameterSet::const_iterator cfgBinningEntry = cfgBinning.begin(); cfgBinningEntry != cfgBinning.end();
0018          ++cfgBinningEntry) {
0019       binning_.emplace_back(new binningEntryType(*cfgBinningEntry));
0020     }
0021   } else {
0022     binning_.emplace_back(new binningEntryType());
0023   }
0024 
0025   for (auto binningEntry = binning_.begin(); binningEntry != binning_.end(); ++binningEntry) {
0026     produces<CorrMETData>((*binningEntry)->binLabel_);
0027   }
0028 }
0029 
0030 PFCandMETcorrInputProducer::~PFCandMETcorrInputProducer() {}
0031 
0032 void PFCandMETcorrInputProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0033   for (auto binningEntry = binning_.begin(); binningEntry != binning_.end(); ++binningEntry) {
0034     (*binningEntry)->binUnclEnergySum_ = CorrMETData();
0035   }
0036 
0037   typedef edm::View<reco::Candidate> CandidateView;
0038   edm::Handle<CandidateView> cands;
0039   evt.getByToken(token_, cands);
0040 
0041   edm::Handle<edm::ValueMap<float>> weights;
0042   if (!weightsToken_.isUninitialized())
0043     evt.getByToken(weightsToken_, weights);
0044 
0045   for (edm::View<reco::Candidate>::const_iterator cand = cands->begin(); cand != cands->end(); ++cand) {
0046     float weight = (!weightsToken_.isUninitialized()) ? (*weights)[cands->ptrAt(cand - cands->begin())] : 1.0;
0047     for (auto binningEntry = binning_.begin(); binningEntry != binning_.end(); ++binningEntry) {
0048       if (!(*binningEntry)->binSelection_ || (*(*binningEntry)->binSelection_)(cand->p4() * weight)) {
0049         (*binningEntry)->binUnclEnergySum_.mex += cand->px() * weight;
0050         (*binningEntry)->binUnclEnergySum_.mey += cand->py() * weight;
0051         (*binningEntry)->binUnclEnergySum_.sumet += cand->et() * weight;
0052       }
0053     }
0054   }
0055 
0056   //--- add momentum sum of PFCandidates not within jets ("unclustered energy") to the event
0057   for (auto binningEntry = binning_.cbegin(); binningEntry != binning_.cend(); ++binningEntry) {
0058     evt.put(std::make_unique<CorrMETData>((*binningEntry)->binUnclEnergySum_), (*binningEntry)->binLabel_);
0059   }
0060 }
0061 
0062 #include "FWCore/Framework/interface/MakerMacros.h"
0063 
0064 DEFINE_FWK_MODULE(PFCandMETcorrInputProducer);