File indexing completed on 2024-04-06 12:26:47
0001 #include "RecoMET/METProducers/interface/ParticleFlowForChargedMETProducer.h"
0002
0003 #include <DataFormats/VertexReco/interface/Vertex.h>
0004 #include <DataFormats/ParticleFlowCandidate/interface/PFCandidate.h>
0005
0006 using namespace edm;
0007 using namespace std;
0008 using namespace reco;
0009
0010 ParticleFlowForChargedMETProducer::ParticleFlowForChargedMETProducer(const edm::ParameterSet& iConfig) {
0011 pfCollectionLabel = iConfig.getParameter<edm::InputTag>("PFCollectionLabel");
0012 pvCollectionLabel = iConfig.getParameter<edm::InputTag>("PVCollectionLabel");
0013
0014 pfCandidatesToken = consumes<PFCandidateCollection>(pfCollectionLabel);
0015 pvCollectionToken = consumes<VertexCollection>(pvCollectionLabel);
0016
0017 dzCut = iConfig.getParameter<double>("dzCut");
0018 neutralEtThreshold = iConfig.getParameter<double>("neutralEtThreshold");
0019
0020 produces<PFCandidateCollection>();
0021 }
0022
0023 void ParticleFlowForChargedMETProducer::produce(Event& iEvent, const EventSetup& iSetup) {
0024
0025 Handle<VertexCollection> pvCollection;
0026 iEvent.getByToken(pvCollectionToken, pvCollection);
0027 VertexCollection::const_iterator vertex = pvCollection->begin();
0028
0029
0030 Handle<PFCandidateCollection> pfCandidates;
0031 iEvent.getByToken(pfCandidatesToken, pfCandidates);
0032
0033
0034 auto chargedPFCandidates = std::make_unique<PFCandidateCollection>();
0035 if (!pvCollection->empty()) {
0036 for (unsigned i = 0; i < pfCandidates->size(); i++) {
0037 const PFCandidate& pfCand = (*pfCandidates)[i];
0038 PFCandidatePtr pfCandPtr(pfCandidates, i);
0039
0040 if (pfCandPtr->trackRef().isNonnull()) {
0041 if (pfCandPtr->trackRef()->dz((*vertex).position()) < dzCut) {
0042 chargedPFCandidates->push_back(pfCand);
0043 chargedPFCandidates->back().setSourceCandidatePtr(pfCandPtr);
0044 }
0045
0046 } else if (neutralEtThreshold > 0 and pfCandPtr->pt() > neutralEtThreshold) {
0047 chargedPFCandidates->push_back(pfCand);
0048 chargedPFCandidates->back().setSourceCandidatePtr(pfCandPtr);
0049 }
0050 }
0051 }
0052
0053 iEvent.put(std::move(chargedPFCandidates));
0054
0055 return;
0056 }
0057
0058 ParticleFlowForChargedMETProducer::~ParticleFlowForChargedMETProducer() {}