Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:31

0001 //
0002 //
0003 
0004 #include "FWCore/Utilities/interface/EDMException.h"
0005 #include "CommonTools/CandUtils/interface/pdgIdUtils.h"
0006 #include "AnalysisDataFormats/TopObjects/interface/StGenEvent.h"
0007 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0008 
0009 StGenEvent::StGenEvent() {}
0010 
0011 StGenEvent::StGenEvent(reco::GenParticleRefProd& parts, reco::GenParticleRefProd& inits) {
0012   parts_ = parts;
0013   initPartons_ = inits;
0014 }
0015 
0016 StGenEvent::~StGenEvent() {}
0017 
0018 const reco::GenParticle* StGenEvent::decayB() const {
0019   const reco::GenParticle* cand = nullptr;
0020   if (singleLepton()) {
0021     const reco::GenParticleCollection& partsColl = *parts_;
0022     const reco::GenParticle& singleLep = *singleLepton();
0023     for (unsigned int i = 0; i < parts_->size(); ++i) {
0024       if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
0025           reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
0026         // ... but it should be the opposite!
0027         cand = &partsColl[i];
0028       }
0029     }
0030   }
0031   return cand;
0032 }
0033 
0034 const reco::GenParticle* StGenEvent::associatedB() const {
0035   const reco::GenParticle* cand = nullptr;
0036   if (singleLepton()) {
0037     const reco::GenParticleCollection& partsColl = *parts_;
0038     const reco::GenParticle& singleLep = *singleLepton();
0039     for (unsigned int i = 0; i < parts_->size(); ++i) {
0040       if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
0041           reco::flavour(singleLep) == reco::flavour(partsColl[i])) {
0042         // ... but it should be the opposite!
0043         cand = &partsColl[i];
0044       }
0045     }
0046   }
0047   return cand;
0048 }
0049 
0050 const reco::GenParticle* StGenEvent::singleLepton() const {
0051   const reco::GenParticle* cand = nullptr;
0052   const reco::GenParticleCollection& partsColl = *parts_;
0053   for (unsigned int i = 0; i < partsColl.size(); ++i) {
0054     if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
0055         std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
0056       cand = &partsColl[i];
0057     }
0058   }
0059   return cand;
0060 }
0061 
0062 const reco::GenParticle* StGenEvent::singleNeutrino() const {
0063   const reco::GenParticle* cand = nullptr;
0064   const reco::GenParticleCollection& partsColl = *parts_;
0065   for (unsigned int i = 0; i < partsColl.size(); ++i) {
0066     if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
0067         std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
0068       cand = &partsColl[i];
0069     }
0070   }
0071   return cand;
0072 }
0073 
0074 const reco::GenParticle* StGenEvent::singleW() const {
0075   const reco::GenParticle* cand = nullptr;
0076   if (singleLepton()) {
0077     const reco::GenParticleCollection& partsColl = *parts_;
0078     const reco::GenParticle& singleLep = *singleLepton();
0079     for (unsigned int i = 0; i < partsColl.size(); ++i) {
0080       if (std::abs(partsColl[i].pdgId()) == TopDecayID::WID &&
0081           reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
0082         // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
0083         cand = &partsColl[i];
0084       }
0085     }
0086   }
0087   return cand;
0088 }
0089 
0090 const reco::GenParticle* StGenEvent::singleTop() const {
0091   const reco::GenParticle* cand = nullptr;
0092   if (singleLepton()) {
0093     const reco::GenParticleCollection& partsColl = *parts_;
0094     const reco::GenParticle& singleLep = *singleLepton();
0095     for (unsigned int i = 0; i < partsColl.size(); ++i) {
0096       if (std::abs(partsColl[i].pdgId()) == TopDecayID::tID &&
0097           reco::flavour(singleLep) != reco::flavour(partsColl[i])) {
0098         cand = &partsColl[i];
0099       }
0100     }
0101   }
0102   return cand;
0103 }