CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
12  pfCollectionLabel = iConfig.getParameter<edm::InputTag>("PFCollectionLabel");
13  pvCollectionLabel = iConfig.getParameter<edm::InputTag>("PVCollectionLabel");
14 
15  pfCandidatesToken = consumes<PFCandidateCollection>(pfCollectionLabel);
16  pvCollectionToken = consumes<VertexCollection>(pvCollectionLabel);
17 
18  dzCut = iConfig.getParameter<double>("dzCut");
19  neutralEtThreshold = iConfig.getParameter<double>("neutralEtThreshold");
20 
21  produces<PFCandidateCollection>();
22 }
23 
24 void ParticleFlowForChargedMETProducer::produce(Event& iEvent, const EventSetup& iSetup)
25 {
26 
27  //Get the PV collection
28  Handle<VertexCollection> pvCollection;
29  iEvent.getByToken(pvCollectionToken, pvCollection);
30  VertexCollection::const_iterator vertex = pvCollection->begin();
31 
32  //Get pfCandidates
34  iEvent.getByToken(pfCandidatesToken, pfCandidates);
35 
36  // the output collection
37  auto_ptr<PFCandidateCollection> chargedPFCandidates( new PFCandidateCollection ) ;
38  if (pvCollection->size()>0) {
39  for( unsigned i=0; i<pfCandidates->size(); i++ ) {
40  const PFCandidate& pfCand = (*pfCandidates)[i];
41  PFCandidatePtr pfCandPtr(pfCandidates, i);
42 
43  if (pfCandPtr->trackRef().isNonnull()) {
44  if (pfCandPtr->trackRef()->dz((*vertex).position()) < dzCut) {
45  chargedPFCandidates->push_back( pfCand );
46  chargedPFCandidates->back().setSourceCandidatePtr( pfCandPtr );
47  }
48 
49  }
50  else if (neutralEtThreshold>0 and
51  pfCandPtr->pt()>neutralEtThreshold) {
52  chargedPFCandidates->push_back( pfCand );
53  chargedPFCandidates->back().setSourceCandidatePtr( pfCandPtr );
54  }
55 
56 
57  }
58  }
59 
60 
61  iEvent.put(chargedPFCandidates);
62 
63  return;
64 }
65 
66 ParticleFlowForChargedMETProducer::~ParticleFlowForChargedMETProducer(){}
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39