CMS 3D CMS Logo

StGenEvent.cc
Go to the documentation of this file.
1 //
2 //
3 
8 
10 
12  parts_ = parts;
13  initPartons_ = inits;
14 }
15 
17 
19  const reco::GenParticle* cand = nullptr;
20  if (singleLepton()) {
21  const reco::GenParticleCollection& partsColl = *parts_;
22  const reco::GenParticle& singleLep = *singleLepton();
23  for (unsigned int i = 0; i < parts_->size(); ++i) {
24  if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
25  reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
26  // ... but it should be the opposite!
27  cand = &partsColl[i];
28  }
29  }
30  }
31  return cand;
32 }
33 
35  const reco::GenParticle* cand = nullptr;
36  if (singleLepton()) {
37  const reco::GenParticleCollection& partsColl = *parts_;
38  const reco::GenParticle& singleLep = *singleLepton();
39  for (unsigned int i = 0; i < parts_->size(); ++i) {
40  if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
41  reco::flavour(singleLep) == reco::flavour(partsColl[i])) {
42  // ... but it should be the opposite!
43  cand = &partsColl[i];
44  }
45  }
46  }
47  return cand;
48 }
49 
51  const reco::GenParticle* cand = nullptr;
52  const reco::GenParticleCollection& partsColl = *parts_;
53  for (unsigned int i = 0; i < partsColl.size(); ++i) {
54  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
55  std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
56  cand = &partsColl[i];
57  }
58  }
59  return cand;
60 }
61 
63  const reco::GenParticle* cand = nullptr;
64  const reco::GenParticleCollection& partsColl = *parts_;
65  for (unsigned int i = 0; i < partsColl.size(); ++i) {
66  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
67  std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
68  cand = &partsColl[i];
69  }
70  }
71  return cand;
72 }
73 
75  const reco::GenParticle* cand = nullptr;
76  if (singleLepton()) {
77  const reco::GenParticleCollection& partsColl = *parts_;
78  const reco::GenParticle& singleLep = *singleLepton();
79  for (unsigned int i = 0; i < partsColl.size(); ++i) {
80  if (std::abs(partsColl[i].pdgId()) == TopDecayID::WID &&
81  reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
82  // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
83  cand = &partsColl[i];
84  }
85  }
86  }
87  return cand;
88 }
89 
91  const reco::GenParticle* cand = nullptr;
92  if (singleLepton()) {
93  const reco::GenParticleCollection& partsColl = *parts_;
94  const reco::GenParticle& singleLep = *singleLepton();
95  for (unsigned int i = 0; i < partsColl.size(); ++i) {
96  if (std::abs(partsColl[i].pdgId()) == TopDecayID::tID &&
97  reco::flavour(singleLep) != reco::flavour(partsColl[i])) {
98  cand = &partsColl[i];
99  }
100  }
101  }
102  return cand;
103 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
const reco::GenParticle * singleNeutrino() const
return single neutrino if available; 0 else
Definition: StGenEvent.cc:62
reco::GenParticleRefProd parts_
reference to the top decay chain (has to be kept in the event!)
Definition: TopGenEvent.h:111
StGenEvent()
empty constructor
Definition: StGenEvent.cc:9
static const int bID
Definition: TopGenEvent.h:13
const reco::GenParticle * singleLepton() const
return single lepton if available; 0 else
Definition: StGenEvent.cc:50
const reco::GenParticle * associatedB() const
return associated b
Definition: StGenEvent.cc:34
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:13
reco::GenParticleRefProd initPartons_
reference to the list of initial partons (has to be kept in the event!)
Definition: TopGenEvent.h:113
static const int tID
Definition: TopGenEvent.h:12
~StGenEvent() override
default destructor
Definition: StGenEvent.cc:16
bool isNeutrino(const Candidate &part)
Definition: pdgIdUtils.h:17
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const reco::GenParticle * singleTop() const
return single Top
Definition: StGenEvent.cc:90
const reco::GenParticle * decayB() const
return decay b
Definition: StGenEvent.cc:18
static const int WID
Definition: TopGenEvent.h:17
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:21
const reco::GenParticle * singleW() const
return single W
Definition: StGenEvent.cc:74