CMS 3D CMS Logo

zMCLeptonDaughters.cc
Go to the documentation of this file.
4 #include <cassert>
5 using namespace std;
6 using namespace reco;
7 
8 pair<const Candidate*, const Candidate*> zMCLeptonDaughters(const Candidate & z, int leptonPdgId) {
9  if(z.numberOfDaughters()<2)
10  throw cms::Exception("RuntimeError") <<
11  "calling helper function reco::zMCLeptonDaughters passing a Z candidate"
12  "with less than 2 daughters (" << z.numberOfDaughters() << ").\n";
13  const Candidate * dau0 = z.daughter(0);
14  const Candidate * dau1 = z.daughter(1);
15  for(size_t i0 = 0; i0 < dau0->numberOfDaughters(); ++i0) {
16  const Candidate * ddau0 = dau0->daughter(i0);
17  if(abs(ddau0->pdgId())==leptonPdgId && ddau0->status()==1) {
18  dau0 = ddau0; break;
19  }
20  }
21  for(size_t i1 = 0; i1 < dau1->numberOfDaughters(); ++i1) {
22  const Candidate * ddau1 = dau1->daughter(i1);
23  if(abs(ddau1->pdgId())==leptonPdgId && ddau1->status()==1) {
24  dau1 = ddau1; break;
25  }
26  }
27  assert(abs(dau0->pdgId())==leptonPdgId && dau0->status()==1);
28  assert(abs(dau1->pdgId())==leptonPdgId && dau1->status()==1);
29  return make_pair(dau0, dau1);
30 }
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual int status() const =0
status word
virtual int pdgId() const =0
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
pair< const Candidate *, const Candidate * > zMCLeptonDaughters(const Candidate &z, int leptonPdgId)
fixed size matrix
virtual size_type numberOfDaughters() const =0
number of daughters