CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GenParticlesHelper.cc
Go to the documentation of this file.
3 
4 using namespace reco;
5 
6 namespace GenParticlesHelper {
7 
8  void findParticles(const reco::GenParticleCollection& sourceParticles,
9  reco::GenParticleRefVector& particleRefs,
10  int pdgId,
11  int status) {
12  unsigned index = 0;
13  for (IG ig = sourceParticles.begin(); ig != sourceParticles.end(); ++ig, ++index) {
14  const GenParticle& gen = *ig;
15 
16  // status has been specified, and this one does not have the correct
17  // status
18  if (status && gen.status() != status)
19  continue;
20 
21  if (std::abs(gen.pdgId()) == pdgId) {
22  GenParticleRef genref(&sourceParticles, index);
23  particleRefs.push_back(genref);
24  }
25  }
26  }
27 
29  reco::GenParticleRefVector& descendents,
30  int status,
31  int pdgId) {
32  const GenParticleRefVector& daughterRefs = base->daughterRefVector();
33 
34  for (IGR idr = daughterRefs.begin(); idr != daughterRefs.end(); ++idr) {
35  if ((*idr)->status() == status && (!pdgId || std::abs((*idr)->pdgId()) == pdgId)) {
36  // cout<<"adding "<<(*idr)<<endl;
37  descendents.push_back(*idr);
38  } else
39  findDescendents(*idr, descendents, status, pdgId);
40  }
41  }
42 
43  void findSisters(const reco::GenParticleRef& baseSister, GenParticleRefVector& sisterRefs) {
44  assert(baseSister->numberOfMothers() > 0);
45 
46  // get first mother
47  const GenParticleRefVector& mothers = baseSister->motherRefVector();
48 
49  // get sisters
50  const GenParticleRefVector allRefs = mothers[0]->daughterRefVector();
51 
53  for (IT id = allRefs.begin(); id != allRefs.end(); ++id) {
54  if (*id == baseSister) {
55  continue; // this is myself
56  } else
57  sisterRefs.push_back(*id);
58  }
59  }
60 
61  bool isDirect(const reco::GenParticleRef& particle) {
62  assert((particle->status() != 0) && (particle->status() < 4));
63  if (particle->status() == 3)
64  return true;
65  else {
66  assert(particle->numberOfMothers() > 0);
67 
68  // get first mother
69  const GenParticleRefVector& mothers = particle->motherRefVector();
70  if (mothers[0]->status() == 3)
71  return true;
72  else
73  return false;
74  }
75  }
76 
77  bool hasAncestor(const reco::GenParticle* particle, int pdgId, int status) {
78  if (particle->pdgId() == pdgId && particle->status() == status)
79  return true;
80 
81  const GenParticleRefVector& mothers = particle->motherRefVector();
82 
83  for (IGR im = mothers.begin(); im != mothers.end(); ++im) {
84  const GenParticle& part = **im;
85  if (hasAncestor(&part, pdgId, status))
86  return true;
87  }
88 
89  return false;
90  }
91 
92  std::ostream& operator<<(std::ostream& out, const reco::GenParticleRef& genRef) {
93  if (!out)
94  return out;
95 
96  out << genRef.key() << " " << genRef->pt();
97 
98  return out;
99  }
100 
101 } // namespace GenParticlesHelper
tuple base
Main Program
Definition: newFWLiteAna.py:92
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
reco::GenParticleRefVector::const_iterator IGR
uint16_t *__restrict__ id
void findParticles(const reco::GenParticleCollection &sourceParticles, reco::GenParticleRefVector &particleRefs, int pdgId, int status)
find all particles of a given pdgId and status
bool isDirect(const reco::GenParticleRef &particle)
check if particle is direct (has status 3 or is a daughter of particle with status 3) ...
list status
Definition: mps_update.py:107
int status() const final
status word
key_type key() const
Accessor for product key.
Definition: Ref.h:250
reco::GenParticleCollection::const_iterator IG
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
assert(be >=bs)
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
int pdgId() const final
PDG identifier.
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
const mothers & motherRefVector() const
references to mothers
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< LinkConnSpec >::const_iterator IT
part
Definition: HCALResponse.h:20
bool hasAncestor(const reco::GenParticle *particle, int pdgId, int status)
does the particle have an ancestor with this pdgId and this status?
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:67
void findSisters(const reco::GenParticleRef &baseSister, reco::GenParticleRefVector &sisterRefs)
find the particles having the same daughter as baseSister
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)