Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 
0003 //____________________________________________________________________________||
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "DataFormats/Candidate/interface/Candidate.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 
0013 #include "DataFormats/METReco/interface/CaloMET.h"
0014 #include "DataFormats/METReco/interface/CorrMETData.h"
0015 
0016 #include "JetMETCorrections/Type1MET/interface/AddCorrectionsToGenericMET.h"
0017 
0018 #include <vector>
0019 
0020 //____________________________________________________________________________||
0021 class CorrectedCaloMETProducer : public edm::stream::EDProducer<> {
0022 public:
0023   explicit CorrectedCaloMETProducer(const edm::ParameterSet& cfg)
0024       : corrector(), token_(consumes<METCollection>(cfg.getParameter<edm::InputTag>("src"))) {
0025     std::vector<edm::InputTag> corrInputTags = cfg.getParameter<std::vector<edm::InputTag> >("srcCorrections");
0026     std::vector<edm::EDGetTokenT<CorrMETData> > corrTokens;
0027     for (std::vector<edm::InputTag>::const_iterator inputTag = corrInputTags.begin(); inputTag != corrInputTags.end();
0028          ++inputTag) {
0029       corrTokens.push_back(consumes<CorrMETData>(*inputTag));
0030     }
0031 
0032     corrector.setCorTokens(corrTokens);
0033 
0034     produces<METCollection>("");
0035   }
0036 
0037   ~CorrectedCaloMETProducer() override {}
0038 
0039 private:
0040   AddCorrectionsToGenericMET corrector;
0041 
0042   typedef std::vector<reco::CaloMET> METCollection;
0043 
0044   edm::EDGetTokenT<METCollection> token_;
0045 
0046   void produce(edm::Event& evt, const edm::EventSetup& es) override {
0047     edm::Handle<METCollection> srcMETCollection;
0048     evt.getByToken(token_, srcMETCollection);
0049 
0050     const reco::CaloMET& srcMET = (*srcMETCollection)[0];
0051 
0052     reco::CaloMET outMET = corrector.getCorrectedCaloMET(srcMET, evt);
0053 
0054     std::unique_ptr<METCollection> product(new METCollection);
0055     product->push_back(outMET);
0056     evt.put(std::move(product));
0057   }
0058 };
0059 
0060 //____________________________________________________________________________||
0061 
0062 DEFINE_FWK_MODULE(CorrectedCaloMETProducer);