00001 // 00002 // $Id: TtGenEvent.cc,v 1.31 2010/10/15 22:44:30 wmtan Exp $ 00003 // 00004 00005 #include "FWCore/Utilities/interface/EDMException.h" 00006 #include "CommonTools/CandUtils/interface/pdgIdUtils.h" 00007 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h" 00008 #include "DataFormats/HepMCCandidate/interface/GenParticle.h" 00009 00010 00012 TtGenEvent::TtGenEvent(reco::GenParticleRefProd& decaySubset, reco::GenParticleRefProd& initSubset) 00013 { 00014 parts_ = decaySubset; 00015 initPartons_= initSubset; 00016 } 00017 00018 WDecay::LeptonType 00019 TtGenEvent::semiLeptonicChannel() const 00020 { 00021 WDecay::LeptonType type=WDecay::kNone; 00022 if( isSemiLeptonic() && singleLepton() ){ 00023 if( std::abs(singleLepton()->pdgId())==TopDecayID::elecID ) type=WDecay::kElec; 00024 if( std::abs(singleLepton()->pdgId())==TopDecayID::muonID ) type=WDecay::kMuon; 00025 if( std::abs(singleLepton()->pdgId())==TopDecayID::tauID ) type=WDecay::kTau; 00026 } 00027 return type; 00028 } 00029 00030 std::pair<WDecay::LeptonType, WDecay::LeptonType> 00031 TtGenEvent::fullLeptonicChannel() const 00032 { 00033 WDecay::LeptonType typeA=WDecay::kNone, typeB=WDecay::kNone; 00034 if( isFullLeptonic() ){ 00035 if( lepton() ){ 00036 if( std::abs(lepton()->pdgId())==TopDecayID::elecID ) typeA=WDecay::kElec; 00037 if( std::abs(lepton()->pdgId())==TopDecayID::muonID ) typeA=WDecay::kMuon; 00038 if( std::abs(lepton()->pdgId())==TopDecayID::tauID ) typeA=WDecay::kTau; 00039 } 00040 if( leptonBar() ){ 00041 if( std::abs(leptonBar()->pdgId())==TopDecayID::elecID ) typeB=WDecay::kElec; 00042 if( std::abs(leptonBar()->pdgId())==TopDecayID::muonID ) typeB=WDecay::kMuon; 00043 if( std::abs(leptonBar()->pdgId())==TopDecayID::tauID ) typeB=WDecay::kTau; 00044 } 00045 } 00046 return ( std::pair<WDecay::LeptonType,WDecay::LeptonType>(typeA, typeB) ); 00047 } 00048 00049 const reco::GenParticle* 00050 TtGenEvent::lepton(bool excludeTauLeptons) const 00051 { 00052 const reco::GenParticle* cand = 0; 00053 const reco::GenParticleCollection& partsColl = *parts_; 00054 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00055 if (reco::isLepton(partsColl[i]) && partsColl[i].mother() && 00056 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00057 if(reco::flavour(partsColl[i])>0){ 00058 cand = &partsColl[i]; 00059 } 00060 } 00061 } 00062 return cand; 00063 } 00064 00065 const reco::GenParticle* 00066 TtGenEvent::leptonBar(bool excludeTauLeptons) const 00067 { 00068 const reco::GenParticle* cand = 0; 00069 const reco::GenParticleCollection& partsColl = *parts_; 00070 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00071 if (reco::isLepton(partsColl[i]) && partsColl[i].mother() && 00072 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00073 if(reco::flavour(partsColl[i])<0){ 00074 cand = &partsColl[i]; 00075 } 00076 } 00077 } 00078 return cand; 00079 } 00080 00081 const reco::GenParticle* 00082 TtGenEvent::singleLepton(bool excludeTauLeptons) const 00083 { 00084 const reco::GenParticle* cand = 0; 00085 if( isSemiLeptonic(excludeTauLeptons) ){ 00086 const reco::GenParticleCollection& partsColl = *parts_; 00087 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00088 if (reco::isLepton(partsColl[i]) && partsColl[i].mother() && 00089 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00090 cand = &partsColl[i]; 00091 } 00092 } 00093 } 00094 return cand; 00095 } 00096 00097 const reco::GenParticle* 00098 TtGenEvent::neutrino(bool excludeTauLeptons) const 00099 { 00100 const reco::GenParticle* cand=0; 00101 const reco::GenParticleCollection & partsColl = *parts_; 00102 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00103 if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() && 00104 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00105 if(reco::flavour(partsColl[i])>0){ 00106 cand = &partsColl[i]; 00107 } 00108 } 00109 } 00110 return cand; 00111 } 00112 00113 const reco::GenParticle* 00114 TtGenEvent::neutrinoBar(bool excludeTauLeptons) const 00115 { 00116 const reco::GenParticle* cand=0; 00117 const reco::GenParticleCollection & partsColl = *parts_; 00118 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00119 if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() && 00120 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00121 if(reco::flavour(partsColl[i])<0){ 00122 cand = &partsColl[i]; 00123 } 00124 } 00125 } 00126 return cand; 00127 } 00128 00129 const reco::GenParticle* 00130 TtGenEvent::singleNeutrino(bool excludeTauLeptons) const 00131 { 00132 const reco::GenParticle* cand=0; 00133 if( isSemiLeptonic(excludeTauLeptons) ) { 00134 const reco::GenParticleCollection & partsColl = *parts_; 00135 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00136 if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() && 00137 std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) { 00138 cand = &partsColl[i]; 00139 } 00140 } 00141 } 00142 return cand; 00143 } 00144 00145 const reco::GenParticle* 00146 TtGenEvent::hadronicDecayQuark(bool invertFlavor) const 00147 { 00148 const reco::GenParticle* cand=0; 00149 // catch W boson and check its daughters for a quark; 00150 // make sure the decay is semi-leptonic first; this 00151 // only makes sense if taus are not excluded from the 00152 // decision 00153 if( singleLepton(false) ){ 00154 for(reco::GenParticleCollection::const_iterator w=parts_->begin(); w!=parts_->end(); ++w){ 00155 if( std::abs( w->pdgId() )==TopDecayID::WID ){ 00156 // make sure that the particle is a W daughter 00157 for(reco::GenParticle::const_iterator wd=w->begin(); wd!=w->end(); ++wd){ 00158 // make sure that the parton is a quark 00159 if( std::abs(wd->pdgId())<TopDecayID::tID ){ 00160 if( invertFlavor?reco::flavour(*wd)<0:reco::flavour(*wd)>0 ){ 00161 cand = dynamic_cast<const reco::GenParticle* > (&(*wd)); 00162 if(cand == 0){ 00163 throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" ); 00164 } 00165 break; 00166 } 00167 } 00168 } 00169 } 00170 } 00171 } 00172 return cand; 00173 } 00174 00175 const reco::GenParticle* 00176 TtGenEvent::hadronicDecayB(bool excludeTauLeptons) const 00177 { 00178 const reco::GenParticle* cand=0; 00179 if( singleLepton(excludeTauLeptons) ){ 00180 const reco::GenParticleCollection& partsColl = *parts_; 00181 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00182 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00183 if (std::abs(partsColl[i].pdgId())==TopDecayID::bID && 00184 reco::flavour(singleLep)==reco::flavour(partsColl[i])) { 00185 cand = &partsColl[i]; 00186 } 00187 } 00188 } 00189 return cand; 00190 } 00191 00192 const reco::GenParticle* 00193 TtGenEvent::hadronicDecayW(bool excludeTauLeptons) const 00194 { 00195 const reco::GenParticle* cand=0; 00196 if( singleLepton(excludeTauLeptons) ){ 00197 const reco::GenParticleCollection& partsColl = *parts_; 00198 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00199 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00200 if (std::abs(partsColl[i].pdgId())==TopDecayID::WID && 00201 reco::flavour(singleLep) != -reco::flavour(partsColl[i])) { 00202 // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign 00203 cand = &partsColl[i]; 00204 } 00205 } 00206 } 00207 return cand; 00208 } 00209 00210 const reco::GenParticle* 00211 TtGenEvent::hadronicDecayTop(bool excludeTauLeptons) const 00212 { 00213 const reco::GenParticle* cand=0; 00214 if( singleLepton(excludeTauLeptons) ){ 00215 const reco::GenParticleCollection& partsColl = *parts_; 00216 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00217 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00218 if (std::abs(partsColl[i].pdgId())==TopDecayID::tID && 00219 reco::flavour(singleLep)==reco::flavour(partsColl[i])) { 00220 cand = &partsColl[i]; 00221 } 00222 } 00223 } 00224 return cand; 00225 } 00226 00227 const reco::GenParticle* 00228 TtGenEvent::leptonicDecayB(bool excludeTauLeptons) const 00229 { 00230 const reco::GenParticle* cand=0; 00231 if( singleLepton(excludeTauLeptons) ){ 00232 const reco::GenParticleCollection& partsColl = *parts_; 00233 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00234 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00235 if (std::abs(partsColl[i].pdgId())==TopDecayID::bID && 00236 reco::flavour(singleLep)!=reco::flavour(partsColl[i])) { 00237 cand = &partsColl[i]; 00238 } 00239 } 00240 } 00241 return cand; 00242 } 00243 00244 const reco::GenParticle* 00245 TtGenEvent::leptonicDecayW(bool excludeTauLeptons) const 00246 { 00247 const reco::GenParticle* cand=0; 00248 if( singleLepton(excludeTauLeptons) ){ 00249 const reco::GenParticleCollection& partsColl = *parts_; 00250 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00251 for (unsigned int i = 0; i < partsColl.size(); ++i) { 00252 if (std::abs(partsColl[i].pdgId())==TopDecayID::WID && 00253 reco::flavour(singleLep) == - reco::flavour(partsColl[i])) { 00254 // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign 00255 cand = &partsColl[i]; 00256 } 00257 } 00258 } 00259 return cand; 00260 } 00261 00262 const reco::GenParticle* 00263 TtGenEvent::leptonicDecayTop(bool excludeTauLeptons) const 00264 { 00265 const reco::GenParticle* cand=0; 00266 if( singleLepton(excludeTauLeptons) ){ 00267 const reco::GenParticleCollection& partsColl = *parts_; 00268 const reco::GenParticle& singleLep = *singleLepton(excludeTauLeptons); 00269 for( unsigned int i = 0; i < partsColl.size(); ++i ){ 00270 if( std::abs(partsColl[i].pdgId())==TopDecayID::tID && 00271 reco::flavour(singleLep)!=reco::flavour(partsColl[i]) ){ 00272 cand = &partsColl[i]; 00273 } 00274 } 00275 } 00276 return cand; 00277 } 00278 00279 std::vector<const reco::GenParticle*> TtGenEvent::leptonicDecayTopRadiation(bool excludeTauLeptons) const{ 00280 if( leptonicDecayTop(excludeTauLeptons) ){ 00281 return (leptonicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID)); 00282 } 00283 std::vector<const reco::GenParticle*> rad; 00284 return (rad); 00285 } 00286 00287 std::vector<const reco::GenParticle*> TtGenEvent::hadronicDecayTopRadiation(bool excludeTauLeptons) const{ 00288 if( hadronicDecayTop(excludeTauLeptons) ){ 00289 return (hadronicDecayTop(excludeTauLeptons)->pdgId()>0 ? radiatedGluons(TopDecayID::tID) : radiatedGluons(-TopDecayID::tID)); 00290 } 00291 std::vector<const reco::GenParticle*> rad; 00292 return (rad); 00293 }