CMS 3D CMS Logo

StGenEvent Class Reference

#include <AnalysisDataFormats/TopObjects/interface/StGenEvent.h>

Inheritance diagram for StGenEvent:

TopGenEvent

List of all members.

Public Member Functions

const reco::GenParticleassociatedB () const
 return associated b
const reco::GenParticledecayB () const
 return decay b
const reco::GenParticlesingleLepton () const
 return single lepton if available; 0 else
const reco::GenParticlesingleNeutrino () const
 return single neutrino if available; 0 else
const reco::GenParticlesingleTop () const
 return single Top
const reco::GenParticlesingleW () const
 return single W
 StGenEvent (reco::GenParticleRefProd &, reco::GenParticleRefProd &)
 default constructor
 StGenEvent ()
 empty constructor
virtual ~StGenEvent ()
 default destructor


Detailed Description

Definition at line 18 of file StGenEvent.h.


Constructor & Destructor Documentation

StGenEvent::StGenEvent (  ) 

empty constructor

Definition at line 10 of file StGenEvent.cc.

00011 {
00012 }

StGenEvent::StGenEvent ( reco::GenParticleRefProd parts,
reco::GenParticleRefProd inits 
)

default constructor

Definition at line 14 of file StGenEvent.cc.

References TopGenEvent::initPartons_, and TopGenEvent::parts_.

00015 {
00016   parts_ = parts;
00017   initPartons_= inits;
00018 }

StGenEvent::~StGenEvent (  )  [virtual]

default destructor

Definition at line 20 of file StGenEvent.cc.

00021 {
00022 }


Member Function Documentation

const reco::GenParticle * StGenEvent::associatedB (  )  const

return associated b

Definition at line 43 of file StGenEvent.cc.

References funct::abs(), TopDecayID::bID, reco::flavour(), i, TopGenEvent::parts_, and singleLepton().

00044 {
00045   const reco::GenParticle* cand=0;
00046   if( singleLepton() ){
00047     const reco::GenParticleCollection & partsColl = *parts_;
00048     const reco::GenParticle & singleLep = *singleLepton();
00049     for (unsigned int i = 0; i < parts_->size(); ++i) {
00050       if (abs(partsColl[i].pdgId())==TopDecayID::bID && 
00051           reco::flavour(singleLep)== reco::flavour(partsColl[i])) { 
00052         // ... but it should be the opposite!
00053         cand = &partsColl[i];
00054       }
00055     }
00056   }
00057   return cand;
00058 }

const reco::GenParticle * StGenEvent::decayB (  )  const

return decay b

Definition at line 25 of file StGenEvent.cc.

References funct::abs(), TopDecayID::bID, reco::flavour(), i, TopGenEvent::parts_, and singleLepton().

00026 {
00027   const reco::GenParticle* cand=0;
00028   if( singleLepton() ){
00029     const reco::GenParticleCollection & partsColl = *parts_;
00030     const reco::GenParticle & singleLep = *singleLepton();
00031     for (unsigned int i = 0; i < parts_->size(); ++i) {
00032       if (abs(partsColl[i].pdgId())==TopDecayID::bID && 
00033           reco::flavour(singleLep)== -reco::flavour(partsColl[i])) { 
00034         // ... but it should be the opposite!
00035         cand = &partsColl[i];
00036       }
00037     }
00038   }
00039   return cand;
00040 }

const reco::GenParticle * StGenEvent::singleLepton (  )  const

return single lepton if available; 0 else

Definition at line 61 of file StGenEvent.cc.

References funct::abs(), i, reco::isLepton(), TopGenEvent::parts_, and TopDecayID::WID.

Referenced by associatedB(), decayB(), singleTop(), and singleW().

00062 {
00063   const reco::GenParticle* cand = 0;
00064   const reco::GenParticleCollection& partsColl = *parts_;
00065   for (unsigned int i = 0; i < partsColl.size(); ++i) {
00066     if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
00067         abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
00068       cand = &partsColl[i];
00069     }
00070   }
00071   return cand;
00072 }

const reco::GenParticle * StGenEvent::singleNeutrino (  )  const

return single neutrino if available; 0 else

Definition at line 75 of file StGenEvent.cc.

References funct::abs(), i, reco::isNeutrino(), TopGenEvent::parts_, and TopDecayID::WID.

00076 {
00077   const reco::GenParticle* cand=0;
00078   const reco::GenParticleCollection & partsColl = *parts_;
00079   for (unsigned int i = 0; i < partsColl.size(); ++i) {
00080     if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
00081         abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
00082       cand = &partsColl[i];
00083     }
00084   }
00085   return cand;
00086 }

const reco::GenParticle * StGenEvent::singleTop (  )  const

return single Top

Definition at line 107 of file StGenEvent.cc.

References funct::abs(), reco::flavour(), i, TopGenEvent::parts_, singleLepton(), and TopDecayID::tID.

00108 {
00109   const reco::GenParticle* cand=0;
00110   if (singleLepton()) {
00111     const reco::GenParticleCollection & partsColl = *parts_;
00112     const reco::GenParticle & singleLep = *singleLepton();
00113     for (unsigned int i = 0; i < partsColl.size(); ++i) {
00114       if (abs(partsColl[i].pdgId())==TopDecayID::tID &&
00115           reco::flavour(singleLep)!=reco::flavour(partsColl[i])) {
00116         cand = &partsColl[i];
00117       }
00118     }
00119   }
00120   return cand;
00121 }

const reco::GenParticle * StGenEvent::singleW (  )  const

return single W

Definition at line 89 of file StGenEvent.cc.

References funct::abs(), reco::flavour(), i, TopGenEvent::parts_, singleLepton(), and TopDecayID::WID.

00090 {
00091   const reco::GenParticle* cand=0;
00092   if (singleLepton()) {
00093     const reco::GenParticleCollection & partsColl = *parts_;
00094     const reco::GenParticle & singleLep = *singleLepton();
00095     for (unsigned int i = 0; i < partsColl.size(); ++i) {
00096       if (abs(partsColl[i].pdgId())==TopDecayID::WID &&
00097           reco::flavour(singleLep) == - reco::flavour(partsColl[i])) { 
00098         // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
00099         cand = &partsColl[i];
00100       }
00101     }
00102   }
00103   return cand;
00104 }


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:32:50 2009 for CMSSW by  doxygen 1.5.4