File indexing completed on 2025-01-31 02:19:56
0001
0002
0003
0004
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 produces<reco::PFClusterMETCollection>().setBranchAlias(iConfig.getParameter<std::string>("alias"));
0036 }
0037
0038
0039 void PFClusterMETProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040 edm::ParameterSetDescription desc;
0041 desc.add<edm::InputTag>("src", edm::InputTag(""));
0042 desc.add<double>("globalThreshold", 0.);
0043 desc.add<std::string>("alias", "");
0044 descriptions.addWithDefaultLabel(desc);
0045 }
0046
0047
0048 void PFClusterMETProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0049 edm::Handle<edm::View<reco::Candidate> > input;
0050 event.getByToken(inputToken_, input);
0051
0052 METAlgo algo;
0053 CommonMETData commonMETdata = algo.run(*input.product(), globalThreshold_);
0054
0055 PFClusterSpecificAlgo pfcluster;
0056 auto pfclustermetcoll = std::make_unique<reco::PFClusterMETCollection>();
0057
0058 pfclustermetcoll->push_back(pfcluster.addInfo(input, commonMETdata));
0059 event.put(std::move(pfclustermetcoll));
0060 }
0061
0062
0063 DEFINE_FWK_MODULE(PFClusterMETProducer);
0064 }
0065
0066