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