CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
TtGenEvent Class Reference

Class derived from the TopGenEvent for ttbar events. More...

#include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"

Inheritance diagram for TtGenEvent:
TopGenEvent

Public Member Functions

bool fromGluonFusion () const
 check if the tops were produced from a pair of gluons More...
 
bool fromQuarkAnnihilation () const
 check if the tops were produced from qqbar More...
 
std::pair< WDecay::LeptonType,
WDecay::LeptonType
fullLeptonicChannel () const
 
const reco::GenParticlehadronicDecayB (bool excludeTauLeptons=false) const
 get b of hadronic decay branch More...
 
const reco::GenParticlehadronicDecayQuark (bool invertFlavor=false) const
 get light quark of hadronic decay branch More...
 
const reco::GenParticlehadronicDecayQuarkBar () const
 get light anti-quark of hadronic decay branch More...
 
const reco::GenParticlehadronicDecayTop (bool excludeTauLeptons=false) const
 get top of hadronic decay branch More...
 
std::vector< const
reco::GenParticle * > 
hadronicDecayTopRadiation (bool excludeTauLeptons=false) const
 gluons as radiated from the hadronicly decaying top quark More...
 
const reco::GenParticlehadronicDecayW (bool excludeTauLeptons=false) const
 get W of hadronic decay branch More...
 
bool isFullHadronic (bool excludeTauLeptons=false) const
 check if the event can be classified as full hadronic More...
 
bool isFullLeptonic (bool excludeTauLeptons=false) const
 check if the event can be classified as full leptonic More...
 
bool isFullLeptonic (WDecay::LeptonType typeA, WDecay::LeptonType typeB) const
 check if the event is full leptonic with the lepton being of typeA or typeB irrelevant of order; all leptons including taus are allowed More...
 
bool isSemiLeptonic (bool excludeTauLeptons=false) const
 check if the event can be classified as semi-laptonic More...
 
bool isSemiLeptonic (WDecay::LeptonType typeA) const
 check if the event is semi-leptonic with the lepton being of typeA; all leptons including taus are allowed More...
 
bool isSemiLeptonic (WDecay::LeptonType typeA, WDecay::LeptonType typeB) const
 check if the event is semi-leptonic with the lepton being of typeA or typeB; all leptons including taus are allowed More...
 
bool isTtBar () const
 check if the event can be classified as ttbar More...
 
const reco::GenParticlelepton (bool excludeTauLeptons=false) const
 get lepton for semi-leptonic or full leptonic decays More...
 
const reco::GenParticleleptonBar (bool excludeTauLeptons=false) const
 get anti-lepton for semi-leptonic or full leptonic decays More...
 
const reco::GenParticleleptonicDecayB (bool excludeTauLeptons=false) const
 get b of leptonic decay branch More...
 
const reco::GenParticleleptonicDecayTop (bool excludeTauLeptons=false) const
 get top of leptonic decay branch More...
 
std::vector< const
reco::GenParticle * > 
leptonicDecayTopRadiation (bool excludeTauLeptons=false) const
 gluons as radiated from the leptonicly decaying top quark More...
 
const reco::GenParticleleptonicDecayW (bool excludeTauLeptons=false) const
 get W of leptonic decay branch More...
 
const reco::GenParticleneutrino (bool excludeTauLeptons=false) const
 get neutrino for semi-leptonic or full leptonic decays More...
 
const reco::GenParticleneutrinoBar (bool excludeTauLeptons=false) const
 get anti-neutrino for semi-leptonic or full leptonic decays More...
 
WDecay::LeptonType semiLeptonicChannel () const
 return decay channel; all leptons including taus are allowed More...
 
const reco::GenParticlesingleLepton (bool excludeTauLeptons=false) const
 return single lepton if available; 0 else More...
 
const reco::GenParticlesingleNeutrino (bool excludeTauLeptons=false) const
 return single neutrino if available; 0 else More...
 
const math::XYZTLorentzVectortopPair () const
 return combined 4-vector of top and topBar More...
 
 TtGenEvent ()
 empty constructor More...
 
 TtGenEvent (reco::GenParticleRefProd &decaySubset, reco::GenParticleRefProd &initSubset)
 default constructor from decaySubset and initSubset More...
 
virtual ~TtGenEvent ()
 default destructor More...
 
- Public Member Functions inherited from TopGenEvent
const reco::GenParticleb () const
 return b quark if available; 0 else More...
 
const reco::GenParticlebBar () const
 return anti-b quark if available; 0 else More...
 
const reco::GenParticlecandidate (int id, unsigned int parentId=0) const
 get candidate with given pdg id if available; 0 else More...
 
const reco::GenParticledaughterQuarkBarOfWMinus () const
 return anti-quark daughter of anti-W boson More...
 
const reco::GenParticledaughterQuarkBarOfWPlus () const
 return anti-quark daughter of W boson More...
 
const reco::GenParticledaughterQuarkOfTop (bool invertCharge=false) const
 return daughter quark of top quark (which can have flavor b, s or d) More...
 
const reco::GenParticledaughterQuarkOfTopBar () const
 return daughter quark of anti-top quark (which can have flavor b, s or d) More...
 
const reco::GenParticledaughterQuarkOfWMinus () const
 return quark daughter of anti-W boson More...
 
const reco::GenParticledaughterQuarkOfWPlus (bool invertQuarkCharge=false, bool invertBosonCharge=false) const
 return quark daughter quark of W boson More...
 
const reco::GenParticleeMinus () const
 return electron if available; 0 else More...
 
const reco::GenParticleePlus () const
 return positron if available; 0 else More...
 
const reco::GenParticleCollectioninitialPartons () const
 return particles of initial partons More...
 
std::vector< const
reco::GenParticle * > 
lightQuarks (bool includingBQuarks=false) const
 return all light quarks or all quarks including b's More...
 
const reco::GenParticlemuMinus () const
 return muon if available; 0 else More...
 
const reco::GenParticlemuPlus () const
 return anti-muon if available; 0 else More...
 
int numberOfBQuarks (bool fromTopQuark=true) const
 return number of b quarks in the decay chain More...
 
int numberOfLeptons (bool fromWBoson=true) const
 return number of leptons in the decay chain More...
 
int numberOfLeptons (WDecay::LeptonType type, bool fromWBoson=true) const
 return number of leptons in the decay chain More...
 
const reco::GenParticleCollectionparticles () const
 return particles of decay chain More...
 
void print () const
 
std::vector< const
reco::GenParticle * > 
radiatedGluons (int pdgId) const
 return radiated gluons from particle with pdgId More...
 
const reco::GenParticletauMinus () const
 return tau if available; 0 else More...
 
const reco::GenParticletauPlus () const
 return anti-tau if available; 0 else More...
 
const reco::GenParticletop () const
 return top if available; 0 else More...
 
const reco::GenParticletopBar () const
 return anti-top if available; 0 else More...
 
 TopGenEvent ()
 empty constructor More...
 
 TopGenEvent (reco::GenParticleRefProd &decaySubset, reco::GenParticleRefProd &iniSubset)
 default constructor More...
 
std::vector< const
reco::GenParticle * > 
topSisters () const
 return number of top anti-top sisters More...
 
const reco::GenParticlewMinus () const
 return W minus if available; 0 else More...
 
const reco::GenParticlewPlus () const
 return W plus if available; 0 else More...
 
virtual ~TopGenEvent ()
 default destructor More...
 

Protected Attributes

math::XYZTLorentzVector topPair_
 combined 4-vector of top and topBar More...
 
- Protected Attributes inherited from TopGenEvent
reco::GenParticleRefProd initPartons_
 reference to the list of initial partons (has to be kept in the event!) More...
 
reco::GenParticleRefProd parts_
 reference to the top decay chain (has to be kept in the event!) More...
 

Private Member Functions

bool isNumberOfLeptons (bool excludeTauLeptons, int nlep) const
 

Detailed Description

Class derived from the TopGenEvent for ttbar events.

The structure holds reference information to the generator particles of the decay chains for each top quark and of the initial partons and provides access and administration. The derived class contains a few additional getters with respect to its base class.

Definition at line 18 of file TtGenEvent.h.

Constructor & Destructor Documentation

TtGenEvent::TtGenEvent ( )
inline

empty constructor

Definition at line 23 of file TtGenEvent.h.

23 {};
TtGenEvent::TtGenEvent ( reco::GenParticleRefProd decaySubset,
reco::GenParticleRefProd initSubset 
)

default constructor from decaySubset and initSubset

Definition at line 10 of file TtGenEvent.cc.

virtual TtGenEvent::~TtGenEvent ( )
inlinevirtual

default destructor

Definition at line 27 of file TtGenEvent.h.

27 {};

Member Function Documentation

bool TtGenEvent::fromGluonFusion ( ) const

check if the tops were produced from a pair of gluons

Definition at line 19 of file TtGenEvent.cc.

bool TtGenEvent::fromQuarkAnnihilation ( ) const

check if the tops were produced from qqbar

Definition at line 29 of file TtGenEvent.cc.

std::pair< WDecay::LeptonType, WDecay::LeptonType > TtGenEvent::fullLeptonicChannel ( ) const

Definition at line 51 of file TtGenEvent.cc.

Referenced by isFullLeptonic().

const reco::GenParticle * TtGenEvent::hadronicDecayB ( bool  excludeTauLeptons = false) const

get b of hadronic decay branch

Definition at line 196 of file TtGenEvent.cc.

Referenced by TtSemiLepEvtPartons::vec().

const reco::GenParticle * TtGenEvent::hadronicDecayQuark ( bool  invertFlavor = false) const

get light quark of hadronic decay branch

Definition at line 166 of file TtGenEvent.cc.

Referenced by hadronicDecayQuarkBar(), and TtSemiLepEvtPartons::vec().

const reco::GenParticle* TtGenEvent::hadronicDecayQuarkBar ( ) const
inline

get light anti-quark of hadronic decay branch

Definition at line 72 of file TtGenEvent.h.

References hadronicDecayQuark().

Referenced by TtSemiLepEvtPartons::vec().

72 {return hadronicDecayQuark(true); };
const reco::GenParticle * hadronicDecayQuark(bool invertFlavor=false) const
get light quark of hadronic decay branch
Definition: TtGenEvent.cc:166
const reco::GenParticle * TtGenEvent::hadronicDecayTop ( bool  excludeTauLeptons = false) const

get top of hadronic decay branch

Definition at line 231 of file TtGenEvent.cc.

std::vector< const reco::GenParticle * > TtGenEvent::hadronicDecayTopRadiation ( bool  excludeTauLeptons = false) const

gluons as radiated from the hadronicly decaying top quark

Definition at line 307 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::hadronicDecayW ( bool  excludeTauLeptons = false) const

get W of hadronic decay branch

Definition at line 213 of file TtGenEvent.cc.

bool TtGenEvent::isFullHadronic ( bool  excludeTauLeptons = false) const
inline

check if the event can be classified as full hadronic

Definition at line 36 of file TtGenEvent.h.

References isNumberOfLeptons(), and isTtBar().

Referenced by TtFullHadEvtPartons::vec().

36 { return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 0) : false;}
bool isNumberOfLeptons(bool excludeTauLeptons, int nlep) const
Definition: TtGenEvent.h:98
volatile std::atomic< bool > shutdown_flag false
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30
bool TtGenEvent::isFullLeptonic ( bool  excludeTauLeptons = false) const
inline

check if the event can be classified as full leptonic

Definition at line 40 of file TtGenEvent.h.

References isNumberOfLeptons(), and isTtBar().

Referenced by TtFullLepEvtPartons::vec().

40 { return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 2) : false;}
bool isNumberOfLeptons(bool excludeTauLeptons, int nlep) const
Definition: TtGenEvent.h:98
volatile std::atomic< bool > shutdown_flag false
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30
bool TtGenEvent::isFullLeptonic ( WDecay::LeptonType  typeA,
WDecay::LeptonType  typeB 
) const
inline

check if the event is full leptonic with the lepton being of typeA or typeB irrelevant of order; all leptons including taus are allowed

Definition at line 102 of file TtGenEvent.h.

References plotBeamSpotDB::first, fullLeptonicChannel(), and edm::second().

103 {
104  return ( (fullLeptonicChannel().first==typeA && fullLeptonicChannel().second==typeB)||
105  (fullLeptonicChannel().first==typeB && fullLeptonicChannel().second==typeA));
106 }
U second(std::pair< T, U > const &p)
std::pair< WDecay::LeptonType, WDecay::LeptonType > fullLeptonicChannel() const
Definition: TtGenEvent.cc:51
bool TtGenEvent::isNumberOfLeptons ( bool  excludeTauLeptons,
int  nlep 
) const
inlineprivate

check whether the number of leptons among the daughters of the W boson is nlep or not; there is an option to exclude taus from the list of leptons to consider

Definition at line 98 of file TtGenEvent.h.

References WDecay::kTau, and TopGenEvent::numberOfLeptons().

Referenced by isFullHadronic(), isFullLeptonic(), and isSemiLeptonic().

98 {return excludeTauLeptons ? (numberOfLeptons()-numberOfLeptons(WDecay::kTau))==nlep : numberOfLeptons()==nlep;}
int numberOfLeptons(bool fromWBoson=true) const
return number of leptons in the decay chain
Definition: TopGenEvent.cc:48
bool TtGenEvent::isSemiLeptonic ( bool  excludeTauLeptons = false) const
inline

check if the event can be classified as semi-laptonic

Definition at line 38 of file TtGenEvent.h.

References isNumberOfLeptons(), and isTtBar().

Referenced by TtSemiLepEvtPartons::vec().

38 { return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 1) : false;}
bool isNumberOfLeptons(bool excludeTauLeptons, int nlep) const
Definition: TtGenEvent.h:98
volatile std::atomic< bool > shutdown_flag false
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30
bool TtGenEvent::isSemiLeptonic ( WDecay::LeptonType  typeA) const
inline

check if the event is semi-leptonic with the lepton being of typeA; all leptons including taus are allowed

Definition at line 45 of file TtGenEvent.h.

References semiLeptonicChannel().

45 { return semiLeptonicChannel()==typeA ? true : false; };
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:39
bool TtGenEvent::isSemiLeptonic ( WDecay::LeptonType  typeA,
WDecay::LeptonType  typeB 
) const
inline

check if the event is semi-leptonic with the lepton being of typeA or typeB; all leptons including taus are allowed

Definition at line 47 of file TtGenEvent.h.

References semiLeptonicChannel().

47 { return (semiLeptonicChannel()==typeA || semiLeptonicChannel()==typeB)? true : false; };
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:39
bool TtGenEvent::isTtBar ( ) const
inline

check if the event can be classified as ttbar

Definition at line 30 of file TtGenEvent.h.

References TopGenEvent::top(), and TopGenEvent::topBar().

Referenced by isFullHadronic(), isFullLeptonic(), and isSemiLeptonic().

30 {return (top() && topBar());}
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:104
const reco::GenParticle * topBar() const
return anti-top if available; 0 else
Definition: TopGenEvent.h:106
const reco::GenParticle * TtGenEvent::lepton ( bool  excludeTauLeptons = false) const

get lepton for semi-leptonic or full leptonic decays

Definition at line 70 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::leptonBar ( bool  excludeTauLeptons = false) const

get anti-lepton for semi-leptonic or full leptonic decays

Definition at line 86 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::leptonicDecayB ( bool  excludeTauLeptons = false) const

get b of leptonic decay branch

Definition at line 248 of file TtGenEvent.cc.

Referenced by TtSemiLepEvtPartons::vec().

const reco::GenParticle * TtGenEvent::leptonicDecayTop ( bool  excludeTauLeptons = false) const

get top of leptonic decay branch

Definition at line 283 of file TtGenEvent.cc.

std::vector< const reco::GenParticle * > TtGenEvent::leptonicDecayTopRadiation ( bool  excludeTauLeptons = false) const

gluons as radiated from the leptonicly decaying top quark

Definition at line 299 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::leptonicDecayW ( bool  excludeTauLeptons = false) const

get W of leptonic decay branch

Definition at line 265 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::neutrino ( bool  excludeTauLeptons = false) const

get neutrino for semi-leptonic or full leptonic decays

Definition at line 118 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::neutrinoBar ( bool  excludeTauLeptons = false) const

get anti-neutrino for semi-leptonic or full leptonic decays

Definition at line 134 of file TtGenEvent.cc.

WDecay::LeptonType TtGenEvent::semiLeptonicChannel ( ) const

return decay channel; all leptons including taus are allowed

Definition at line 39 of file TtGenEvent.cc.

Referenced by isSemiLeptonic().

const reco::GenParticle * TtGenEvent::singleLepton ( bool  excludeTauLeptons = false) const

return single lepton if available; 0 else

Definition at line 102 of file TtGenEvent.cc.

const reco::GenParticle * TtGenEvent::singleNeutrino ( bool  excludeTauLeptons = false) const

return single neutrino if available; 0 else

Definition at line 150 of file TtGenEvent.cc.

const math::XYZTLorentzVector* TtGenEvent::topPair ( ) const
inline

return combined 4-vector of top and topBar

Definition at line 87 of file TtGenEvent.h.

87 { return isTtBar() ? &topPair_ : 0; };
math::XYZTLorentzVector topPair_
combined 4-vector of top and topBar
Definition: TtGenEvent.h:87
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30

Member Data Documentation

math::XYZTLorentzVector TtGenEvent::topPair_
protected

combined 4-vector of top and topBar

Definition at line 87 of file TtGenEvent.h.