CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/AnalysisDataFormats/TopObjects/src/TopGenEvent.cc

Go to the documentation of this file.
00001 #include "FWCore/Utilities/interface/EDMException.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 #include "CommonTools/CandUtils/interface/pdgIdUtils.h"
00005 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00006 #include "AnalysisDataFormats/TopObjects/interface/TopGenEvent.h"
00007 
00009 TopGenEvent::TopGenEvent(reco::GenParticleRefProd& decaySubset, reco::GenParticleRefProd& initSubset)
00010 {
00011   parts_ = decaySubset; 
00012   initPartons_= initSubset;
00013 }
00014 
00015 const reco::GenParticle*
00016 TopGenEvent::candidate(int id, unsigned int parentId) const
00017 {
00018   const reco::GenParticle* cand=0;
00019   const reco::GenParticleCollection & partsColl = *parts_;
00020   for( unsigned int i = 0; i < partsColl.size(); ++i ) {
00021     if( partsColl[i].pdgId()==id ){
00022       if(parentId==0?true:partsColl[i].mother()&&std::abs(partsColl[i].mother()->pdgId())==(int)parentId){
00023         cand = &partsColl[i];
00024       }
00025     }
00026   }  
00027   return cand;
00028 }
00029 
00030 void
00031 TopGenEvent::print() const 
00032 {
00033   edm::LogVerbatim log("TopGenEvent");
00034   log << "\n"
00035       << "--------------------------------------\n"
00036       << "- Dump TopGenEvent Content           -\n"
00037       << "--------------------------------------\n";
00038   for (reco::GenParticleCollection::const_iterator part = parts_->begin(); 
00039        part<parts_->end(); ++part) {
00040     log << "pdgId:"  << std::setw(5)  << part->pdgId()     << ", "
00041         << "mass:"   << std::setw(11) << part->p4().mass() << ", "
00042         << "energy:" << std::setw(11) << part->energy()    << ", " 
00043         << "pt:"     << std::setw(11) << part->pt()        << "\n"; 
00044   }
00045 }
00046 
00047 int
00048 TopGenEvent::numberOfLeptons(bool fromWBoson) const
00049 {
00050   int lep=0;
00051   const reco::GenParticleCollection& partsColl = *parts_;
00052   for(unsigned int i = 0; i < partsColl.size(); ++i) {
00053     if(reco::isLepton(partsColl[i])) {
00054       if(fromWBoson){
00055         if(partsColl[i].mother() &&  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID){
00056           ++lep;
00057         }
00058       }
00059       else{
00060         ++lep;
00061       }
00062     }
00063   }  
00064   return lep;
00065 }
00066 
00067 int
00068 TopGenEvent::numberOfLeptons(WDecay::LeptonType typeRestriction, bool fromWBoson) const
00069 {
00070   int leptonType=-1;
00071   switch(typeRestriction){
00072     // resolve whether or not there is
00073     // any restriction in lepton types
00074   case WDecay::kElec: 
00075     leptonType=TopDecayID::elecID;
00076     break;
00077   case WDecay::kMuon: 
00078     leptonType=TopDecayID::muonID;
00079     break;
00080   case WDecay::kTau: 
00081     leptonType=TopDecayID::tauID;
00082     break;
00083   case WDecay::kNone:
00084     break;
00085   }
00086   int lep=0;
00087   const reco::GenParticleCollection & partsColl = *parts_;
00088   for(unsigned int i = 0; i < partsColl.size(); ++i) {
00089     if(fromWBoson){
00090       // restrict to particles originating from the W boson
00091       if( !(partsColl[i].mother() &&  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::WID) ){
00092         continue;
00093       }
00094     }
00095     if(leptonType>0){
00096       // in case of lepton type restriction
00097       if( std::abs(partsColl[i].pdgId())==leptonType ){
00098         ++lep;
00099       }
00100     }
00101     else{
00102       // take any lepton type into account else
00103       if( reco::isLepton(partsColl[i]) ){
00104         ++lep;
00105       }
00106     }
00107   }  
00108   return lep;
00109 }
00110 
00111 int
00112 TopGenEvent::numberOfBQuarks(bool fromTopQuark) const
00113 {
00114   int bq=0;
00115   const reco::GenParticleCollection & partsColl = *parts_;
00116   for (unsigned int i = 0; i < partsColl.size(); ++i) {
00117    //depend if radiation qqbar are included or not
00118     if(std::abs(partsColl[i].pdgId())==TopDecayID::bID){
00119       if(fromTopQuark){
00120         if(partsColl[i].mother() &&  std::abs(partsColl[i].mother()->pdgId())==TopDecayID::tID){
00121           ++bq;
00122         }
00123       }
00124       else{
00125         ++bq;
00126       }
00127     }
00128   }  
00129   return bq;
00130 }
00131 
00132 std::vector<const reco::GenParticle*> 
00133 TopGenEvent::topSisters() const
00134 {
00135   std::vector<const reco::GenParticle*> sisters;
00136   for(reco::GenParticleCollection::const_iterator part = parts_->begin(); part<parts_->end(); ++part){
00137     if( part->numberOfMothers()==0 && std::abs(part->pdgId())!= TopDecayID::tID){
00138       // choose top sister which do not have a 
00139       // mother and are whether top nor anti-top 
00140       if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == 0){
00141         throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00142       }
00143       sisters.push_back( part->clone() );
00144     }
00145   }  
00146   return sisters;
00147 }
00148 
00149 const reco::GenParticle*
00150 TopGenEvent::daughterQuarkOfTop(bool invertCharge) const
00151 {
00152   const reco::GenParticle* cand=0;
00153   for(reco::GenParticleCollection::const_iterator top = parts_->begin(); top<parts_->end(); ++top){
00154     if( top->pdgId()==(invertCharge?-TopDecayID::tID:TopDecayID::tID) ){
00155       for(reco::GenParticle::const_iterator quark = top->begin(); quark<top->end(); ++quark){
00156         if( std::abs(quark->pdgId())<= TopDecayID::bID ){
00157           cand = dynamic_cast<const reco::GenParticle* > (&(*quark));
00158           if(cand == 0){
00159             throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00160           }
00161           break;
00162         }
00163       }
00164     }
00165   }
00166   return cand;
00167 }
00168 
00169 const reco::GenParticle* 
00170 TopGenEvent::daughterQuarkOfWPlus(bool invertQuarkCharge, bool invertBosonCharge) const 
00171 {
00172   const reco::GenParticle* cand=0;
00173   const reco::GenParticleCollection & partsColl = *parts_;
00174   for (unsigned int i = 0; i < partsColl.size(); ++i) {
00175     if(partsColl[i].mother() && partsColl[i].mother()->pdgId()==(invertBosonCharge?-TopDecayID::WID:TopDecayID::WID) &&
00176        std::abs(partsColl[i].pdgId())<=TopDecayID::bID && (invertQuarkCharge?reco::flavour(partsColl[i])<0:reco::flavour(partsColl[i])>0)){
00177       cand = &partsColl[i];
00178     }
00179   }
00180   return cand;
00181 }
00182 
00183 std::vector<const reco::GenParticle*> 
00184 TopGenEvent::lightQuarks(bool includingBQuarks) const 
00185 {
00186   std::vector<const reco::GenParticle*> lightQuarks;
00187   for (reco::GenParticleCollection::const_iterator part = parts_->begin(); part < parts_->end(); ++part) {
00188     if( (includingBQuarks && std::abs(part->pdgId())==TopDecayID::bID) || std::abs(part->pdgId())<TopDecayID::bID ) {
00189       if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == 0){
00190         throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00191       }
00192       lightQuarks.push_back( part->clone() );
00193     }
00194   }  
00195   return lightQuarks;
00196 }
00197 
00198 std::vector<const reco::GenParticle*> 
00199 TopGenEvent::radiatedGluons(int pdgId) const{
00200   std::vector<const reco::GenParticle*> rads;
00201   for (reco::GenParticleCollection::const_iterator part = parts_->begin(); part < parts_->end(); ++part) {
00202     if ( part->mother() && part->mother()->pdgId()==pdgId ){
00203       if(part->pdgId()==TopDecayID::glueID){
00204         if( dynamic_cast<const reco::GenParticle*>( &(*part) ) == 0){
00205           throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00206         }
00207       }
00208       rads.push_back( part->clone() );
00209     }
00210   }  
00211   return rads;
00212 }