CMS 3D CMS Logo

TtGenEvent.h
Go to the documentation of this file.
1 #ifndef TopObjects_TtGenEvent_h
2 #define TopObjects_TtGenEvent_h
3 
6 
18 class TtGenEvent : public TopGenEvent {
19 public:
25  ~TtGenEvent() override{};
26 
28  bool isTtBar() const { return (top() && topBar()); }
30  bool fromGluonFusion() const;
32  bool fromQuarkAnnihilation() const;
34  bool isFullHadronic(bool excludeTauLeptons = false) const {
35  return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 0) : false;
36  }
38  bool isSemiLeptonic(bool excludeTauLeptons = false) const {
39  return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 1) : false;
40  }
42  bool isFullLeptonic(bool excludeTauLeptons = false) const {
43  return isTtBar() ? isNumberOfLeptons(excludeTauLeptons, 2) : false;
44  }
45 
49  bool isSemiLeptonic(WDecay::LeptonType typeA) const { return semiLeptonicChannel() == typeA ? true : false; };
52  return (semiLeptonicChannel() == typeA || semiLeptonicChannel() == typeB) ? true : false;
53  };
54  // return decay channel (as a std::pair of LeptonType's); all leptons including taus are allowed
55  std::pair<WDecay::LeptonType, WDecay::LeptonType> fullLeptonicChannel() const;
57  bool isFullLeptonic(WDecay::LeptonType typeA, WDecay::LeptonType typeB) const;
58 
60  const reco::GenParticle* singleLepton(bool excludeTauLeptons = false) const;
62  const reco::GenParticle* singleNeutrino(bool excludeTauLeptons = false) const;
64  const reco::GenParticle* leptonicDecayW(bool excludeTauLeptons = false) const;
66  const reco::GenParticle* leptonicDecayB(bool excludeTauLeptons = false) const;
68  const reco::GenParticle* leptonicDecayTop(bool excludeTauLeptons = false) const;
70  const reco::GenParticle* hadronicDecayW(bool excludeTauLeptons = false) const;
72  const reco::GenParticle* hadronicDecayB(bool excludeTauLeptons = false) const;
74  const reco::GenParticle* hadronicDecayTop(bool excludeTauLeptons = false) const;
76  const reco::GenParticle* hadronicDecayQuark(bool invertFlavor = false) const;
80  std::vector<const reco::GenParticle*> leptonicDecayTopRadiation(bool excludeTauLeptons = false) const;
82  std::vector<const reco::GenParticle*> hadronicDecayTopRadiation(bool excludeTauLeptons = false) const;
84  const reco::GenParticle* lepton(bool excludeTauLeptons = false) const;
86  const reco::GenParticle* leptonBar(bool excludeTauLeptons = false) const;
88  const reco::GenParticle* neutrino(bool excludeTauLeptons = false) const;
90  const reco::GenParticle* neutrinoBar(bool excludeTauLeptons = false) const;
91 
93  const math::XYZTLorentzVector* topPair() const { return isTtBar() ? &topPair_ : nullptr; };
94 
95 protected:
98 
99 private:
102  bool isNumberOfLeptons(bool excludeTauLeptons, int nlep) const {
103  return excludeTauLeptons ? (numberOfLeptons() - numberOfLeptons(WDecay::kTau)) == nlep : numberOfLeptons() == nlep;
104  }
105 };
106 
108  return ((fullLeptonicChannel().first == typeA && fullLeptonicChannel().second == typeB) ||
109  (fullLeptonicChannel().first == typeB && fullLeptonicChannel().second == typeA));
110 }
111 
112 #endif
bool isNumberOfLeptons(bool excludeTauLeptons, int nlep) const
Definition: TtGenEvent.h:102
std::vector< const reco::GenParticle * > leptonicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the leptonicly decaying top quark
Definition: TtGenEvent.cc:273
const reco::GenParticle * hadronicDecayTop(bool excludeTauLeptons=false) const
get top of hadronic decay branch
Definition: TtGenEvent.cc:212
bool fromQuarkAnnihilation() const
check if the tops were produced from qqbar
Definition: TtGenEvent.cc:25
bool isFullHadronic(bool excludeTauLeptons=false) const
check if the event can be classified as full hadronic
Definition: TtGenEvent.h:34
bool isSemiLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as semi-laptonic
Definition: TtGenEvent.h:38
const reco::GenParticle * leptonicDecayB(bool excludeTauLeptons=false) const
get b of leptonic decay branch
Definition: TtGenEvent.cc:227
int numberOfLeptons(bool fromWBoson=true) const
return number of leptons in the decay chain
Definition: TopGenEvent.cc:41
math::XYZTLorentzVector topPair_
combined 4-vector of top and topBar
Definition: TtGenEvent.h:93
const reco::GenParticle * hadronicDecayB(bool excludeTauLeptons=false) const
get b of hadronic decay branch
Definition: TtGenEvent.cc:181
const reco::GenParticle * topBar() const
return anti-top if available; 0 else
Definition: TopGenEvent.h:103
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
std::vector< const reco::GenParticle * > hadronicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the hadronicly decaying top quark
Definition: TtGenEvent.cc:282
const reco::GenParticle * hadronicDecayQuarkBar() const
get light anti-quark of hadronic decay branch
Definition: TtGenEvent.h:78
std::pair< WDecay::LeptonType, WDecay::LeptonType > fullLeptonicChannel() const
Definition: TtGenEvent.cc:46
TtGenEvent()
empty constructor
Definition: TtGenEvent.h:21
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:101
const reco::GenParticle * leptonBar(bool excludeTauLeptons=false) const
get anti-lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:83
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:28
const reco::GenParticle * singleLepton(bool excludeTauLeptons=false) const
return single lepton if available; 0 else
Definition: TtGenEvent.cc:97
const reco::GenParticle * hadronicDecayQuark(bool invertFlavor=false) const
get light quark of hadronic decay branch
Definition: TtGenEvent.cc:153
const reco::GenParticle * hadronicDecayW(bool excludeTauLeptons=false) const
get W of hadronic decay branch
Definition: TtGenEvent.cc:196
const reco::GenParticle * neutrino(bool excludeTauLeptons=false) const
get neutrino for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:111
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 ta...
Definition: TtGenEvent.h:51
Base class to hold information for reduced top generator information.
Definition: TopGenEvent.h:40
bool isSemiLeptonic(WDecay::LeptonType typeA) const
check if the event is semi-leptonic with the lepton being of typeA; all leptons including taus are al...
Definition: TtGenEvent.h:49
bool isFullLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as full leptonic
Definition: TtGenEvent.h:42
const reco::GenParticle * lepton(bool excludeTauLeptons=false) const
get lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:69
WDecay::LeptonType semiLeptonicChannel() const
return decay channel; all leptons including taus are allowed
Definition: TtGenEvent.cc:33
const reco::GenParticle * neutrinoBar(bool excludeTauLeptons=false) const
get anti-neutrino for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:125
const reco::GenParticle * singleNeutrino(bool excludeTauLeptons=false) const
return single neutrino if available; 0 else
Definition: TtGenEvent.cc:139
const reco::GenParticle * leptonicDecayW(bool excludeTauLeptons=false) const
get W of leptonic decay branch
Definition: TtGenEvent.cc:242
bool fromGluonFusion() const
check if the tops were produced from a pair of gluons
Definition: TtGenEvent.cc:17
const reco::GenParticle * leptonicDecayTop(bool excludeTauLeptons=false) const
get top of leptonic decay branch
Definition: TtGenEvent.cc:258
const math::XYZTLorentzVector * topPair() const
return combined 4-vector of top and topBar
Definition: TtGenEvent.h:93
~TtGenEvent() override
default destructor
Definition: TtGenEvent.h:25