CMS 3D CMS Logo

ParticleFlowForChargedMETProducer.cc
Go to the documentation of this file.
2 
5 
6 using namespace edm;
7 using namespace std;
8 using namespace reco;
9 
10 ParticleFlowForChargedMETProducer::ParticleFlowForChargedMETProducer(const edm::ParameterSet& iConfig) {
11  pfCollectionLabel = iConfig.getParameter<edm::InputTag>("PFCollectionLabel");
12  pvCollectionLabel = iConfig.getParameter<edm::InputTag>("PVCollectionLabel");
13 
14  pfCandidatesToken = consumes<PFCandidateCollection>(pfCollectionLabel);
15  pvCollectionToken = consumes<VertexCollection>(pvCollectionLabel);
16 
17  dzCut = iConfig.getParameter<double>("dzCut");
18  neutralEtThreshold = iConfig.getParameter<double>("neutralEtThreshold");
19 
20  produces<PFCandidateCollection>();
21 }
22 
23 void ParticleFlowForChargedMETProducer::produce(Event& iEvent, const EventSetup& iSetup) {
24  //Get the PV collection
26  iEvent.getByToken(pvCollectionToken, pvCollection);
27  VertexCollection::const_iterator vertex = pvCollection->begin();
28 
29  //Get pfCandidates
31  iEvent.getByToken(pfCandidatesToken, pfCandidates);
32 
33  // the output collection
34  auto chargedPFCandidates = std::make_unique<PFCandidateCollection>();
35  if (!pvCollection->empty()) {
36  for (unsigned i = 0; i < pfCandidates->size(); i++) {
37  const PFCandidate& pfCand = (*pfCandidates)[i];
38  PFCandidatePtr pfCandPtr(pfCandidates, i);
39 
40  if (pfCandPtr->trackRef().isNonnull()) {
41  if (pfCandPtr->trackRef()->dz((*vertex).position()) < dzCut) {
42  chargedPFCandidates->push_back(pfCand);
43  chargedPFCandidates->back().setSourceCandidatePtr(pfCandPtr);
44  }
45 
46  } else if (neutralEtThreshold > 0 and pfCandPtr->pt() > neutralEtThreshold) {
47  chargedPFCandidates->push_back(pfCand);
48  chargedPFCandidates->back().setSourceCandidatePtr(pfCandPtr);
49  }
50  }
51  }
52 
54 
55  return;
56 }
57 
58 ParticleFlowForChargedMETProducer::~ParticleFlowForChargedMETProducer() {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double pt() const final
transverse momentum
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
int iEvent
Definition: GenABIO.cc:224
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
HLT enums.
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:430
def move(src, dest)
Definition: eostools.py:511