Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    METProducers
0004 // Class:      GenMETProducer
0005 //
0006 //
0007 
0008 //____________________________________________________________________________||
0009 #include "RecoMET/METProducers/interface/GenMETProducer.h"
0010 
0011 #include "FWCore/Framework/interface/ConsumesCollector.h"
0012 
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 
0015 #include "DataFormats/Common/interface/Handle.h"
0016 
0017 #include "DataFormats/METReco/interface/METFwd.h"
0018 #include "DataFormats/METReco/interface/GenMETFwd.h"
0019 #include "DataFormats/METReco/interface/CommonMETData.h"
0020 
0021 #include "RecoMET/METAlgorithms/interface/GenSpecificAlgo.h"
0022 
0023 #include <cstring>
0024 
0025 //____________________________________________________________________________||
0026 namespace cms {
0027 
0028   //____________________________________________________________________________||
0029   GenMETProducer::GenMETProducer(const edm::ParameterSet& iConfig)
0030       : inputToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("src"))),
0031         globalThreshold_(iConfig.getParameter<double>("globalThreshold")),
0032         onlyFiducial_(iConfig.getParameter<bool>("onlyFiducialParticles")),
0033         applyFiducialThresholdForFractions_(iConfig.getParameter<bool>("applyFiducialThresholdForFractions")),
0034         usePt_(iConfig.getParameter<bool>("usePt")) {
0035     std::string alias = iConfig.exists("alias") ? iConfig.getParameter<std::string>("alias") : "";
0036     produces<reco::GenMETCollection>().setBranchAlias(alias);
0037   }
0038 
0039   //____________________________________________________________________________||
0040   void GenMETProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0041     edm::Handle<edm::View<reco::Candidate> > input;
0042     event.getByToken(inputToken_, input);
0043 
0044     CommonMETData commonMETdata;
0045 
0046     GenSpecificAlgo gen;
0047     auto genmetcoll = std::make_unique<reco::GenMETCollection>();
0048     genmetcoll->push_back(gen.addInfo(
0049         input, &commonMETdata, globalThreshold_, onlyFiducial_, applyFiducialThresholdForFractions_, usePt_));
0050     event.put(std::move(genmetcoll));
0051   }
0052 
0053   //____________________________________________________________________________||
0054   DEFINE_FWK_MODULE(GenMETProducer);
0055 }  // namespace cms
0056 
0057 //____________________________________________________________________________||