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