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