Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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