CommonTools
CandUtils
src
zMCLeptonDaughters.cc
Go to the documentation of this file.
1
#include "
CommonTools/CandUtils/interface/zMCLeptonDaughters.h
"
2
#include "
DataFormats/Candidate/interface/Candidate.h
"
3
#include "
FWCore/Utilities/interface/EDMException.h
"
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"
) <<
"calling helper function reco::zMCLeptonDaughters passing a Z candidate"
11
"with less than 2 daughters ("
12
<<
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;
19
break
;
20
}
21
}
22
for
(
size_t
i1
= 0;
i1
< dau1->
numberOfDaughters
(); ++
i1
) {
23
const
Candidate
* ddau1 = dau1->
daughter
(
i1
);
24
if
(
abs
(ddau1->
pdgId
()) ==
leptonPdgId
&& ddau1->
status
() == 1) {
25
dau1 = ddau1;
26
break
;
27
}
28
}
29
assert
(
abs
(dau0->
pdgId
()) ==
leptonPdgId
&& dau0->
status
() == 1);
30
assert
(
abs
(dau1->
pdgId
()) ==
leptonPdgId
&& dau1->
status
() == 1);
31
return
make_pair(dau0, dau1);
32
}
reco::Candidate::daughter
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode)
testProducerWithPsetDescEmpty_cfi.i1
i1
Definition:
testProducerWithPsetDescEmpty_cfi.py:45
cms::cuda::assert
assert(be >=bs)
reco::Candidate::status
virtual int status() const =0
status word
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
EDMException.h
DDAxes::z
reco::Candidate::numberOfDaughters
virtual size_type numberOfDaughters() const =0
number of daughters
zMCLeptonDaughters.h
reco::Candidate::pdgId
virtual int pdgId() const =0
PDG identifier.
reco::Candidate
Definition:
Candidate.h:27
std
Definition:
JetResolutionObject.h:76
jets_cff.leptonPdgId
leptonPdgId
Definition:
jets_cff.py:164
cms::Exception
Definition:
Exception.h:70
Candidate.h
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
zMCLeptonDaughters
pair< const Candidate *, const Candidate * > zMCLeptonDaughters(const Candidate &z, int leptonPdgId)
Definition:
zMCLeptonDaughters.cc:8
Generated for CMSSW Reference Manual by
1.8.16