CMS 3D CMS Logo

TopGenEvent.cc
Go to the documentation of this file.
3 
7 
10 {
13 }
14 
15 const reco::GenParticle*
16 TopGenEvent::candidate(int id, unsigned int parentId) const
17 {
18  const reco::GenParticle* cand=nullptr;
19  const reco::GenParticleCollection & partsColl = *parts_;
20  for( unsigned int i = 0; i < partsColl.size(); ++i ) {
21  if( partsColl[i].pdgId()==id ){
22  if(parentId==0?true:partsColl[i].mother()&&std::abs(partsColl[i].mother()->pdgId())==(int)parentId){
23  cand = &partsColl[i];
24  }
25  }
26  }
27  return cand;
28 }
29 
30 void
32 {
33  edm::LogVerbatim log("TopGenEvent");
34  log << "\n"
35  << "--------------------------------------\n"
36  << "- Dump TopGenEvent Content -\n"
37  << "--------------------------------------\n";
38  for (reco::GenParticleCollection::const_iterator part = parts_->begin();
39  part<parts_->end(); ++part) {
40  log << "pdgId:" << std::setw(5) << part->pdgId() << ", "
41  << "mass:" << std::setw(11) << part->p4().mass() << ", "
42  << "energy:" << std::setw(11) << part->energy() << ", "
43  << "pt:" << std::setw(11) << part->pt() << "\n";
44  }
45 }
46 
47 int
48 TopGenEvent::numberOfLeptons(bool fromWBoson) const
49 {
50  int lep=0;
51  const reco::GenParticleCollection& partsColl = *parts_;
52  for(unsigned int i = 0; i < partsColl.size(); ++i) {
53  if(reco::isLepton(partsColl[i])) {
54  if(fromWBoson){
55  if(partsColl[i].mother() && std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID){
56  ++lep;
57  }
58  }
59  else{
60  ++lep;
61  }
62  }
63  }
64  return lep;
65 }
66 
67 int
68 TopGenEvent::numberOfLeptons(WDecay::LeptonType typeRestriction, bool fromWBoson) const
69 {
70  int leptonType=-1;
71  switch(typeRestriction){
72  // resolve whether or not there is
73  // any restriction in lepton types
74  case WDecay::kElec:
75  leptonType=TopDecayID::elecID;
76  break;
77  case WDecay::kMuon:
78  leptonType=TopDecayID::muonID;
79  break;
80  case WDecay::kTau:
81  leptonType=TopDecayID::tauID;
82  break;
83  case WDecay::kNone:
84  break;
85  }
86  int lep=0;
87  const reco::GenParticleCollection & partsColl = *parts_;
88  for(unsigned int i = 0; i < partsColl.size(); ++i) {
89  if(fromWBoson){
90  // restrict to particles originating from the W boson
91  if( !(partsColl[i].mother() && std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) ){
92  continue;
93  }
94  }
95  if(leptonType>0){
96  // in case of lepton type restriction
97  if( std::abs(partsColl[i].pdgId())==leptonType ){
98  ++lep;
99  }
100  }
101  else{
102  // take any lepton type into account else
103  if( reco::isLepton(partsColl[i]) ){
104  ++lep;
105  }
106  }
107  }
108  return lep;
109 }
110 
111 int
112 TopGenEvent::numberOfBQuarks(bool fromTopQuark) const
113 {
114  int bq=0;
115  const reco::GenParticleCollection & partsColl = *parts_;
116  for (unsigned int i = 0; i < partsColl.size(); ++i) {
117  //depend if radiation qqbar are included or not
118  if(std::abs(partsColl[i].pdgId())==TopDecayID::bID){
119  if(fromTopQuark){
120  if(partsColl[i].mother() && std::abs(partsColl[i].mother()->pdgId())==TopDecayID::tID){
121  ++bq;
122  }
123  }
124  else{
125  ++bq;
126  }
127  }
128  }
129  return bq;
130 }
131 
132 std::vector<const reco::GenParticle*>
134 {
135  std::vector<const reco::GenParticle*> sisters;
136  for(reco::GenParticleCollection::const_iterator part = parts_->begin(); part<parts_->end(); ++part){
137  if( part->numberOfMothers()==0 && std::abs(part->pdgId())!= TopDecayID::tID){
138  // choose top sister which do not have a
139  // mother and are whether top nor anti-top
140  if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == nullptr){
141  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
142  }
143  sisters.push_back( part->clone() );
144  }
145  }
146  return sisters;
147 }
148 
149 const reco::GenParticle*
150 TopGenEvent::daughterQuarkOfTop(bool invertCharge) const
151 {
152  const reco::GenParticle* cand=nullptr;
153  for(reco::GenParticleCollection::const_iterator top = parts_->begin(); top<parts_->end(); ++top){
154  if( top->pdgId()==(invertCharge?-TopDecayID::tID:TopDecayID::tID) ){
155  for(reco::GenParticle::const_iterator quark = top->begin(); quark<top->end(); ++quark){
156  if( std::abs(quark->pdgId())<= TopDecayID::bID ){
157  cand = dynamic_cast<const reco::GenParticle* > (&(*quark));
158  if(cand == nullptr){
159  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
160  }
161  break;
162  }
163  }
164  }
165  }
166  return cand;
167 }
168 
169 const reco::GenParticle*
170 TopGenEvent::daughterQuarkOfWPlus(bool invertQuarkCharge, bool invertBosonCharge) const
171 {
172  const reco::GenParticle* cand=nullptr;
173  const reco::GenParticleCollection & partsColl = *parts_;
174  for (unsigned int i = 0; i < partsColl.size(); ++i) {
175  if(partsColl[i].mother() && partsColl[i].mother()->pdgId()==(invertBosonCharge?-TopDecayID::WID:TopDecayID::WID) &&
176  std::abs(partsColl[i].pdgId())<=TopDecayID::bID && (invertQuarkCharge?reco::flavour(partsColl[i])<0:reco::flavour(partsColl[i])>0)){
177  cand = &partsColl[i];
178  }
179  }
180  return cand;
181 }
182 
183 std::vector<const reco::GenParticle*>
184 TopGenEvent::lightQuarks(bool includingBQuarks) const
185 {
186  std::vector<const reco::GenParticle*> lightQuarks;
187  for (reco::GenParticleCollection::const_iterator part = parts_->begin(); part < parts_->end(); ++part) {
188  if( (includingBQuarks && std::abs(part->pdgId())==TopDecayID::bID) || std::abs(part->pdgId())<TopDecayID::bID ) {
189  if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == nullptr){
190  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
191  }
192  lightQuarks.push_back( part->clone() );
193  }
194  }
195  return lightQuarks;
196 }
197 
198 std::vector<const reco::GenParticle*>
200  std::vector<const reco::GenParticle*> rads;
201  for (reco::GenParticleCollection::const_iterator part = parts_->begin(); part < parts_->end(); ++part) {
202  if ( part->mother() && part->mother()->pdgId()==pdgId ){
203  if(part->pdgId()==TopDecayID::glueID){
204  if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == nullptr){
205  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
206  }
207  }
208  rads.push_back( part->clone() );
209  }
210  }
211  return rads;
212 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int pdgId() const final
PDG identifier.
reco::GenParticleRefProd parts_
reference to the top decay chain (has to be kept in the event!)
Definition: TopGenEvent.h:115
static const int bID
Definition: TopGenEvent.h:14
static const int glueID
Definition: TopGenEvent.h:15
std::vector< const reco::GenParticle * > topSisters() const
return number of top anti-top sisters
Definition: TopGenEvent.cc:133
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
TopGenEvent()
empty constructor
Definition: TopGenEvent.h:46
const reco::GenParticle * daughterQuarkOfWPlus(bool invertQuarkCharge=false, bool invertBosonCharge=false) const
return quark daughter quark of W boson
Definition: TopGenEvent.cc:170
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:104
const reco::GenParticle * daughterQuarkOfTop(bool invertCharge=false) const
return daughter quark of top quark (which can have flavor b, s or d)
Definition: TopGenEvent.cc:150
reco::GenParticleRefProd initPartons_
reference to the list of initial partons (has to be kept in the event!)
Definition: TopGenEvent.h:117
static const int tID
Definition: TopGenEvent.h:13
static const int tauID
Definition: TopGenEvent.h:21
leptonType
LEPTON
Definition: autophobj.py:48
const_iterator end() const
last daughter const_iterator
Definition: Candidate.h:146
std::vector< const reco::GenParticle * > lightQuarks(bool includingBQuarks=false) const
return all light quarks or all quarks including b&#39;s
Definition: TopGenEvent.cc:184
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const reco::GenParticle * candidate(int id, unsigned int parentId=0) const
get candidate with given pdg id if available; 0 else
Definition: TopGenEvent.cc:16
int numberOfBQuarks(bool fromTopQuark=true) const
return number of b quarks in the decay chain
Definition: TopGenEvent.cc:112
static const int muonID
Definition: TopGenEvent.h:20
part
Definition: HCALResponse.h:20
std::vector< const reco::GenParticle * > radiatedGluons(int pdgId) const
return radiated gluons from particle with pdgId
Definition: TopGenEvent.cc:199
static const int WID
Definition: TopGenEvent.h:18
static const int elecID
Definition: TopGenEvent.h:19
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:144
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
void print() const
Definition: TopGenEvent.cc:31
int numberOfLeptons(bool fromWBoson=true) const
return number of leptons in the decay chain
Definition: TopGenEvent.cc:48