CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtGenEvent.cc
Go to the documentation of this file.
1 //
2 // $Id: TtGenEvent.cc,v 1.30 2010/01/05 22:33:10 rwolf Exp $
3 //
4 
9 
10 
13 {
16 }
17 
20 {
22  if( isSemiLeptonic() && singleLepton() ){
26  }
27  return type;
28 }
29 
30 std::pair<WDecay::LeptonType, WDecay::LeptonType>
32 {
34  if( isFullLeptonic() ){
35  if( lepton() ){
39  }
40  if( leptonBar() ){
44  }
45  }
46  return ( std::pair<WDecay::LeptonType,WDecay::LeptonType>(typeA, typeB) );
47 }
48 
49 const reco::GenParticle*
50 TtGenEvent::lepton(bool excludeTauLeptons) const
51 {
52  const reco::GenParticle* cand = 0;
53  const reco::GenParticleCollection& partsColl = *parts_;
54  for (unsigned int i = 0; i < partsColl.size(); ++i) {
55  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
56  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
57  if(reco::flavour(partsColl[i])>0){
58  cand = &partsColl[i];
59  }
60  }
61  }
62  return cand;
63 }
64 
65 const reco::GenParticle*
66 TtGenEvent::leptonBar(bool excludeTauLeptons) const
67 {
68  const reco::GenParticle* cand = 0;
69  const reco::GenParticleCollection& partsColl = *parts_;
70  for (unsigned int i = 0; i < partsColl.size(); ++i) {
71  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
72  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
73  if(reco::flavour(partsColl[i])<0){
74  cand = &partsColl[i];
75  }
76  }
77  }
78  return cand;
79 }
80 
81 const reco::GenParticle*
82 TtGenEvent::singleLepton(bool excludeTauLeptons) const
83 {
84  const reco::GenParticle* cand = 0;
85  if( isSemiLeptonic(excludeTauLeptons) ){
86  const reco::GenParticleCollection& partsColl = *parts_;
87  for (unsigned int i = 0; i < partsColl.size(); ++i) {
88  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
89  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
90  cand = &partsColl[i];
91  }
92  }
93  }
94  return cand;
95 }
96 
97 const reco::GenParticle*
98 TtGenEvent::neutrino(bool excludeTauLeptons) const
99 {
100  const reco::GenParticle* cand=0;
101  const reco::GenParticleCollection & partsColl = *parts_;
102  for (unsigned int i = 0; i < partsColl.size(); ++i) {
103  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
104  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
105  if(reco::flavour(partsColl[i])>0){
106  cand = &partsColl[i];
107  }
108  }
109  }
110  return cand;
111 }
112 
113 const reco::GenParticle*
114 TtGenEvent::neutrinoBar(bool excludeTauLeptons) const
115 {
116  const reco::GenParticle* cand=0;
117  const reco::GenParticleCollection & partsColl = *parts_;
118  for (unsigned int i = 0; i < partsColl.size(); ++i) {
119  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
120  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
121  if(reco::flavour(partsColl[i])<0){
122  cand = &partsColl[i];
123  }
124  }
125  }
126  return cand;
127 }
128 
129 const reco::GenParticle*
130 TtGenEvent::singleNeutrino(bool excludeTauLeptons) const
131 {
132  const reco::GenParticle* cand=0;
133  if( isSemiLeptonic(excludeTauLeptons) ) {
134  const reco::GenParticleCollection & partsColl = *parts_;
135  for (unsigned int i = 0; i < partsColl.size(); ++i) {
136  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
137  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
138  cand = &partsColl[i];
139  }
140  }
141  }
142  return cand;
143 }
144 
145 const reco::GenParticle*
146 TtGenEvent::hadronicDecayQuark(bool invertFlavor) const
147 {
148  const reco::GenParticle* cand=0;
149  // catch W boson and check its daughters for a quark;
150  // make sure the decay is semi-leptonic first; this
151  // only makes sense if taus are not excluded from the
152  // decision
153  if( singleLepton(false) ){
154  for(reco::GenParticleCollection::const_iterator w=parts_->begin(); w!=parts_->end(); ++w){
155  if( std::abs( w->pdgId() )==TopDecayID::WID ){
156  // make sure that the particle is a W daughter
157  for(reco::GenParticle::const_iterator wd=w->begin(); wd!=w->end(); ++wd){
158  // make sure that the parton is a quark
159  if( std::abs(wd->pdgId())<TopDecayID::tID ){
160  if( invertFlavor?reco::flavour(*wd)<0:reco::flavour(*wd)>0 ){
161  cand = dynamic_cast<const reco::GenParticle* > (&(*wd));
162  if(cand == 0){
163  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
164  }
165  break;
166  }
167  }
168  }
169  }
170  }
171  }
172  return cand;
173 }
174 
175 const reco::GenParticle*
176 TtGenEvent::hadronicDecayB(bool excludeTauLeptons) const
177 {
178  const reco::GenParticle* cand=0;
179  if( singleLepton(excludeTauLeptons) ){
180  const reco::GenParticleCollection& partsColl = *parts_;
181  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
182  for (unsigned int i = 0; i < partsColl.size(); ++i) {
183  if (std::abs(partsColl[i].pdgId())==TopDecayID::bID &&
184  reco::flavour(singleLep)==reco::flavour(partsColl[i])) {
185  cand = &partsColl[i];
186  }
187  }
188  }
189  return cand;
190 }
191 
192 const reco::GenParticle*
193 TtGenEvent::hadronicDecayW(bool excludeTauLeptons) const
194 {
195  const reco::GenParticle* cand=0;
196  if( singleLepton(excludeTauLeptons) ){
197  const reco::GenParticleCollection& partsColl = *parts_;
198  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
199  for (unsigned int i = 0; i < partsColl.size(); ++i) {
200  if (std::abs(partsColl[i].pdgId())==TopDecayID::WID &&
201  reco::flavour(singleLep) != -reco::flavour(partsColl[i])) {
202  // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
203  cand = &partsColl[i];
204  }
205  }
206  }
207  return cand;
208 }
209 
210 const reco::GenParticle*
211 TtGenEvent::hadronicDecayTop(bool excludeTauLeptons) const
212 {
213  const reco::GenParticle* cand=0;
214  if( singleLepton(excludeTauLeptons) ){
215  const reco::GenParticleCollection& partsColl = *parts_;
216  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
217  for (unsigned int i = 0; i < partsColl.size(); ++i) {
218  if (std::abs(partsColl[i].pdgId())==TopDecayID::tID &&
219  reco::flavour(singleLep)==reco::flavour(partsColl[i])) {
220  cand = &partsColl[i];
221  }
222  }
223  }
224  return cand;
225 }
226 
227 const reco::GenParticle*
228 TtGenEvent::leptonicDecayB(bool excludeTauLeptons) const
229 {
230  const reco::GenParticle* cand=0;
231  if( singleLepton(excludeTauLeptons) ){
232  const reco::GenParticleCollection& partsColl = *parts_;
233  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
234  for (unsigned int i = 0; i < partsColl.size(); ++i) {
235  if (std::abs(partsColl[i].pdgId())==TopDecayID::bID &&
236  reco::flavour(singleLep)!=reco::flavour(partsColl[i])) {
237  cand = &partsColl[i];
238  }
239  }
240  }
241  return cand;
242 }
243 
244 const reco::GenParticle*
245 TtGenEvent::leptonicDecayW(bool excludeTauLeptons) const
246 {
247  const reco::GenParticle* cand=0;
248  if( singleLepton(excludeTauLeptons) ){
249  const reco::GenParticleCollection& partsColl = *parts_;
250  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
251  for (unsigned int i = 0; i < partsColl.size(); ++i) {
252  if (std::abs(partsColl[i].pdgId())==TopDecayID::WID &&
253  reco::flavour(singleLep) == - reco::flavour(partsColl[i])) {
254  // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
255  cand = &partsColl[i];
256  }
257  }
258  }
259  return cand;
260 }
261 
262 const reco::GenParticle*
263 TtGenEvent::leptonicDecayTop(bool excludeTauLeptons) const
264 {
265  const reco::GenParticle* cand=0;
266  if( singleLepton(excludeTauLeptons) ){
267  const reco::GenParticleCollection& partsColl = *parts_;
268  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
269  for( unsigned int i = 0; i < partsColl.size(); ++i ){
270  if( std::abs(partsColl[i].pdgId())==TopDecayID::tID &&
271  reco::flavour(singleLep)!=reco::flavour(partsColl[i]) ){
272  cand = &partsColl[i];
273  }
274  }
275  }
276  return cand;
277 }
278 
279 std::vector<const reco::GenParticle*> TtGenEvent::leptonicDecayTopRadiation(bool excludeTauLeptons) const{
280  if( leptonicDecayTop(excludeTauLeptons) ){
281  return (leptonicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID));
282  }
283  std::vector<const reco::GenParticle*> rad;
284  return (rad);
285 }
286 
287 std::vector<const reco::GenParticle*> TtGenEvent::hadronicDecayTopRadiation(bool excludeTauLeptons) const{
288  if( hadronicDecayTop(excludeTauLeptons) ){
289  return (hadronicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID));
290  }
291  std::vector<const reco::GenParticle*> rad;
292  return (rad);
293 }
type
Definition: HCALResponse.h:22
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int i
Definition: DBlmapReader.cc:9
const reco::GenParticle * leptonicDecayTop(bool excludeTauLeptons=false) const
get top of leptonic decay branch
Definition: TtGenEvent.cc:263
const reco::GenParticle * hadronicDecayB(bool excludeTauLeptons=false) const
get b of hadronic decay branch
Definition: TtGenEvent.cc:176
reco::GenParticleRefProd parts_
reference to the top decay chain (has to be kept in the event!)
Definition: TopGenEvent.h:115
bool isSemiLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as semi-laptonic
Definition: TtGenEvent.h:34
static const int bID
Definition: TopGenEvent.h:14
const reco::GenParticle * lepton(bool excludeTauLeptons=false) const
get lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:50
#define abs(x)
Definition: mlp_lapack.h:159
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
std::vector< const reco::GenParticle * > leptonicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the leptonicly decaying top quark
Definition: TtGenEvent.cc:279
const reco::GenParticle * hadronicDecayW(bool excludeTauLeptons=false) const
get W of hadronic decay branch
Definition: TtGenEvent.cc:193
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
bool isFullLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as full leptonic
Definition: TtGenEvent.h:36
bool isNeutrino(const Candidate &part)
Definition: pdgIdUtils.h:25
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:19
virtual const_iterator end() const =0
last daughter const_iterator
std::vector< const reco::GenParticle * > hadronicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the hadronicly decaying top quark
Definition: TtGenEvent.cc:287
TtGenEvent()
empty constructor
Definition: TtGenEvent.h:23
const reco::GenParticle * hadronicDecayTop(bool excludeTauLeptons=false) const
get top of hadronic decay branch
Definition: TtGenEvent.cc:211
const reco::GenParticle * hadronicDecayQuark(bool invertFlavor=false) const
get light quark of hadronic decay branch
Definition: TtGenEvent.cc:146
std::pair< WDecay::LeptonType, WDecay::LeptonType > fullLeptonicChannel() const
Definition: TtGenEvent.cc:31
const reco::GenParticle * singleNeutrino(bool excludeTauLeptons=false) const
return single neutrino if available; 0 else
Definition: TtGenEvent.cc:130
static const int muonID
Definition: TopGenEvent.h:20
const reco::GenParticle * leptonicDecayB(bool excludeTauLeptons=false) const
get b of leptonic decay branch
Definition: TtGenEvent.cc:228
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
virtual const_iterator begin() const =0
first daughter const_iterator
const reco::GenParticle * neutrino(bool excludeTauLeptons=false) const
get neutrino for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:98
static const int elecID
Definition: TopGenEvent.h:19
const reco::GenParticle * neutrinoBar(bool excludeTauLeptons=false) const
get anti-neutrino for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:114
const reco::GenParticle * leptonicDecayW(bool excludeTauLeptons=false) const
get W of leptonic decay branch
Definition: TtGenEvent.cc:245
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
const reco::GenParticle * leptonBar(bool excludeTauLeptons=false) const
get anti-lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:66
const reco::GenParticle * singleLepton(bool excludeTauLeptons=false) const
return single lepton if available; 0 else
Definition: TtGenEvent.cc:82