CMS 3D CMS Logo

TtGenEvent.cc
Go to the documentation of this file.
1 //
2 //
3 
8 
11 {
14  if(top() && topBar())
16 }
17 
18 bool
20 {
21  const reco::GenParticleCollection& initPartsColl = *initPartons_;
22  if(initPartsColl.size()==2)
23  if(initPartsColl[0].pdgId()==21 && initPartsColl[1].pdgId()==21)
24  return true;
25  return false;
26 }
27 
28 bool
30 {
31  const reco::GenParticleCollection& initPartsColl = *initPartons_;
32  if(initPartsColl.size()==2)
33  if(std::abs(initPartsColl[0].pdgId())<TopDecayID::tID && initPartsColl[0].pdgId()==-initPartsColl[1].pdgId())
34  return true;
35  return false;
36 }
37 
40 {
42  if( isSemiLeptonic() && singleLepton() ){
46  }
47  return type;
48 }
49 
50 std::pair<WDecay::LeptonType, WDecay::LeptonType>
52 {
54  if( isFullLeptonic() ){
55  if( lepton() ){
59  }
60  if( leptonBar() ){
64  }
65  }
66  return ( std::pair<WDecay::LeptonType,WDecay::LeptonType>(typeA, typeB) );
67 }
68 
69 const reco::GenParticle*
70 TtGenEvent::lepton(bool excludeTauLeptons) const
71 {
72  const reco::GenParticle* cand = 0;
73  const reco::GenParticleCollection& partsColl = *parts_;
74  for (unsigned int i = 0; i < partsColl.size(); ++i) {
75  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
76  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
77  if(reco::flavour(partsColl[i])>0){
78  cand = &partsColl[i];
79  }
80  }
81  }
82  return cand;
83 }
84 
85 const reco::GenParticle*
86 TtGenEvent::leptonBar(bool excludeTauLeptons) const
87 {
88  const reco::GenParticle* cand = 0;
89  const reco::GenParticleCollection& partsColl = *parts_;
90  for (unsigned int i = 0; i < partsColl.size(); ++i) {
91  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
92  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
93  if(reco::flavour(partsColl[i])<0){
94  cand = &partsColl[i];
95  }
96  }
97  }
98  return cand;
99 }
100 
101 const reco::GenParticle*
102 TtGenEvent::singleLepton(bool excludeTauLeptons) const
103 {
104  const reco::GenParticle* cand = 0;
105  if( isSemiLeptonic(excludeTauLeptons) ){
106  const reco::GenParticleCollection& partsColl = *parts_;
107  for (unsigned int i = 0; i < partsColl.size(); ++i) {
108  if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
109  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
110  cand = &partsColl[i];
111  }
112  }
113  }
114  return cand;
115 }
116 
117 const reco::GenParticle*
118 TtGenEvent::neutrino(bool excludeTauLeptons) const
119 {
120  const reco::GenParticle* cand=0;
121  const reco::GenParticleCollection & partsColl = *parts_;
122  for (unsigned int i = 0; i < partsColl.size(); ++i) {
123  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
124  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
125  if(reco::flavour(partsColl[i])>0){
126  cand = &partsColl[i];
127  }
128  }
129  }
130  return cand;
131 }
132 
133 const reco::GenParticle*
134 TtGenEvent::neutrinoBar(bool excludeTauLeptons) const
135 {
136  const reco::GenParticle* cand=0;
137  const reco::GenParticleCollection & partsColl = *parts_;
138  for (unsigned int i = 0; i < partsColl.size(); ++i) {
139  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
140  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
141  if(reco::flavour(partsColl[i])<0){
142  cand = &partsColl[i];
143  }
144  }
145  }
146  return cand;
147 }
148 
149 const reco::GenParticle*
150 TtGenEvent::singleNeutrino(bool excludeTauLeptons) const
151 {
152  const reco::GenParticle* cand=0;
153  if( isSemiLeptonic(excludeTauLeptons) ) {
154  const reco::GenParticleCollection & partsColl = *parts_;
155  for (unsigned int i = 0; i < partsColl.size(); ++i) {
156  if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
157  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) {
158  cand = &partsColl[i];
159  }
160  }
161  }
162  return cand;
163 }
164 
165 const reco::GenParticle*
166 TtGenEvent::hadronicDecayQuark(bool invertFlavor) const
167 {
168  const reco::GenParticle* cand=0;
169  // catch W boson and check its daughters for a quark;
170  // make sure the decay is semi-leptonic first; this
171  // only makes sense if taus are not excluded from the
172  // decision
173  if( singleLepton(false) ){
174  for(reco::GenParticleCollection::const_iterator w=parts_->begin(); w!=parts_->end(); ++w){
175  if( std::abs( w->pdgId() )==TopDecayID::WID ){
176  // make sure that the particle is a W daughter
177  for(reco::GenParticle::const_iterator wd=w->begin(); wd!=w->end(); ++wd){
178  // make sure that the parton is a quark
179  if( std::abs(wd->pdgId())<TopDecayID::tID ){
180  if( invertFlavor?reco::flavour(*wd)<0:reco::flavour(*wd)>0 ){
181  cand = dynamic_cast<const reco::GenParticle* > (&(*wd));
182  if(cand == 0){
183  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
184  }
185  break;
186  }
187  }
188  }
189  }
190  }
191  }
192  return cand;
193 }
194 
195 const reco::GenParticle*
196 TtGenEvent::hadronicDecayB(bool excludeTauLeptons) const
197 {
198  const reco::GenParticle* cand=0;
199  if( singleLepton(excludeTauLeptons) ){
200  const reco::GenParticleCollection& partsColl = *parts_;
201  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
202  for (unsigned int i = 0; i < partsColl.size(); ++i) {
203  if (std::abs(partsColl[i].pdgId())==TopDecayID::bID &&
204  reco::flavour(singleLep)==reco::flavour(partsColl[i])) {
205  cand = &partsColl[i];
206  }
207  }
208  }
209  return cand;
210 }
211 
212 const reco::GenParticle*
213 TtGenEvent::hadronicDecayW(bool excludeTauLeptons) const
214 {
215  const reco::GenParticle* cand=0;
216  if( singleLepton(excludeTauLeptons) ){
217  const reco::GenParticleCollection& partsColl = *parts_;
218  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
219  for (unsigned int i = 0; i < partsColl.size(); ++i) {
220  if (std::abs(partsColl[i].pdgId())==TopDecayID::WID &&
221  reco::flavour(singleLep) != -reco::flavour(partsColl[i])) {
222  // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
223  cand = &partsColl[i];
224  }
225  }
226  }
227  return cand;
228 }
229 
230 const reco::GenParticle*
231 TtGenEvent::hadronicDecayTop(bool excludeTauLeptons) const
232 {
233  const reco::GenParticle* cand=0;
234  if( singleLepton(excludeTauLeptons) ){
235  const reco::GenParticleCollection& partsColl = *parts_;
236  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
237  for (unsigned int i = 0; i < partsColl.size(); ++i) {
238  if (std::abs(partsColl[i].pdgId())==TopDecayID::tID &&
239  reco::flavour(singleLep)==reco::flavour(partsColl[i])) {
240  cand = &partsColl[i];
241  }
242  }
243  }
244  return cand;
245 }
246 
247 const reco::GenParticle*
248 TtGenEvent::leptonicDecayB(bool excludeTauLeptons) const
249 {
250  const reco::GenParticle* cand=0;
251  if( singleLepton(excludeTauLeptons) ){
252  const reco::GenParticleCollection& partsColl = *parts_;
253  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
254  for (unsigned int i = 0; i < partsColl.size(); ++i) {
255  if (std::abs(partsColl[i].pdgId())==TopDecayID::bID &&
256  reco::flavour(singleLep)!=reco::flavour(partsColl[i])) {
257  cand = &partsColl[i];
258  }
259  }
260  }
261  return cand;
262 }
263 
264 const reco::GenParticle*
265 TtGenEvent::leptonicDecayW(bool excludeTauLeptons) const
266 {
267  const reco::GenParticle* cand=0;
268  if( singleLepton(excludeTauLeptons) ){
269  const reco::GenParticleCollection& partsColl = *parts_;
270  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
271  for (unsigned int i = 0; i < partsColl.size(); ++i) {
272  if (std::abs(partsColl[i].pdgId())==TopDecayID::WID &&
273  reco::flavour(singleLep) == - reco::flavour(partsColl[i])) {
274  // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
275  cand = &partsColl[i];
276  }
277  }
278  }
279  return cand;
280 }
281 
282 const reco::GenParticle*
283 TtGenEvent::leptonicDecayTop(bool excludeTauLeptons) const
284 {
285  const reco::GenParticle* cand=0;
286  if( singleLepton(excludeTauLeptons) ){
287  const reco::GenParticleCollection& partsColl = *parts_;
288  const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons);
289  for( unsigned int i = 0; i < partsColl.size(); ++i ){
290  if( std::abs(partsColl[i].pdgId())==TopDecayID::tID &&
291  reco::flavour(singleLep)!=reco::flavour(partsColl[i]) ){
292  cand = &partsColl[i];
293  }
294  }
295  }
296  return cand;
297 }
298 
299 std::vector<const reco::GenParticle*> TtGenEvent::leptonicDecayTopRadiation(bool excludeTauLeptons) const{
300  if( leptonicDecayTop(excludeTauLeptons) ){
301  return (leptonicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID));
302  }
303  std::vector<const reco::GenParticle*> rad;
304  return (rad);
305 }
306 
307 std::vector<const reco::GenParticle*> TtGenEvent::hadronicDecayTopRadiation(bool excludeTauLeptons) const{
308  if( hadronicDecayTop(excludeTauLeptons) ){
309  return (hadronicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID));
310  }
311  std::vector<const reco::GenParticle*> rad;
312  return (rad);
313 }
type
Definition: HCALResponse.h:21
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:283
const reco::GenParticle * hadronicDecayB(bool excludeTauLeptons=false) const
get b of hadronic decay branch
Definition: TtGenEvent.cc:196
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:38
static const int bID
Definition: TopGenEvent.h:14
const double w
Definition: UKUtility.cc:23
bool fromQuarkAnnihilation() const
check if the tops were produced from qqbar
Definition: TtGenEvent.cc:29
const reco::GenParticle * lepton(bool excludeTauLeptons=false) const
get lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:70
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:19
math::XYZTLorentzVector topPair_
combined 4-vector of top and topBar
Definition: TtGenEvent.h:87
std::vector< const reco::GenParticle * > leptonicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the leptonicly decaying top quark
Definition: TtGenEvent.cc:299
const reco::GenParticle * hadronicDecayW(bool excludeTauLeptons=false) const
get W of hadronic decay branch
Definition: TtGenEvent.cc:213
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:104
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
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
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:40
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:39
bool fromGluonFusion() const
check if the tops were produced from a pair of gluons
Definition: TtGenEvent.cc:19
double p4[4]
Definition: TauolaWrapper.h:92
std::vector< const reco::GenParticle * > hadronicDecayTopRadiation(bool excludeTauLeptons=false) const
gluons as radiated from the hadronicly decaying top quark
Definition: TtGenEvent.cc:307
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TtGenEvent()
empty constructor
Definition: TtGenEvent.h:23
const reco::GenParticle * hadronicDecayTop(bool excludeTauLeptons=false) const
get top of hadronic decay branch
Definition: TtGenEvent.cc:231
const reco::GenParticle * hadronicDecayQuark(bool invertFlavor=false) const
get light quark of hadronic decay branch
Definition: TtGenEvent.cc:166
std::pair< WDecay::LeptonType, WDecay::LeptonType > fullLeptonicChannel() const
Definition: TtGenEvent.cc:51
const reco::GenParticle * singleNeutrino(bool excludeTauLeptons=false) const
return single neutrino if available; 0 else
Definition: TtGenEvent.cc:150
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:248
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
const reco::GenParticle * neutrino(bool excludeTauLeptons=false) const
get neutrino for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:118
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:134
const reco::GenParticle * leptonicDecayW(bool excludeTauLeptons=false) const
get W of leptonic decay branch
Definition: TtGenEvent.cc:265
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
const reco::GenParticle * topBar() const
return anti-top if available; 0 else
Definition: TopGenEvent.h:106
const reco::GenParticle * leptonBar(bool excludeTauLeptons=false) const
get anti-lepton for semi-leptonic or full leptonic decays
Definition: TtGenEvent.cc:86
const reco::GenParticle * singleLepton(bool excludeTauLeptons=false) const
return single lepton if available; 0 else
Definition: TtGenEvent.cc:102