Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103
//
//

#include "FWCore/Utilities/interface/EDMException.h"
#include "CommonTools/CandUtils/interface/pdgIdUtils.h"
#include "AnalysisDataFormats/TopObjects/interface/StGenEvent.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"

StGenEvent::StGenEvent() {}

StGenEvent::StGenEvent(reco::GenParticleRefProd& parts, reco::GenParticleRefProd& inits) {
  parts_ = parts;
  initPartons_ = inits;
}

StGenEvent::~StGenEvent() {}

const reco::GenParticle* StGenEvent::decayB() const {
  const reco::GenParticle* cand = nullptr;
  if (singleLepton()) {
    const reco::GenParticleCollection& partsColl = *parts_;
    const reco::GenParticle& singleLep = *singleLepton();
    for (unsigned int i = 0; i < parts_->size(); ++i) {
      if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
          reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
        // ... but it should be the opposite!
        cand = &partsColl[i];
      }
    }
  }
  return cand;
}

const reco::GenParticle* StGenEvent::associatedB() const {
  const reco::GenParticle* cand = nullptr;
  if (singleLepton()) {
    const reco::GenParticleCollection& partsColl = *parts_;
    const reco::GenParticle& singleLep = *singleLepton();
    for (unsigned int i = 0; i < parts_->size(); ++i) {
      if (std::abs(partsColl[i].pdgId()) == TopDecayID::bID &&
          reco::flavour(singleLep) == reco::flavour(partsColl[i])) {
        // ... but it should be the opposite!
        cand = &partsColl[i];
      }
    }
  }
  return cand;
}

const reco::GenParticle* StGenEvent::singleLepton() const {
  const reco::GenParticle* cand = nullptr;
  const reco::GenParticleCollection& partsColl = *parts_;
  for (unsigned int i = 0; i < partsColl.size(); ++i) {
    if (reco::isLepton(partsColl[i]) && partsColl[i].mother() &&
        std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
      cand = &partsColl[i];
    }
  }
  return cand;
}

const reco::GenParticle* StGenEvent::singleNeutrino() const {
  const reco::GenParticle* cand = nullptr;
  const reco::GenParticleCollection& partsColl = *parts_;
  for (unsigned int i = 0; i < partsColl.size(); ++i) {
    if (reco::isNeutrino(partsColl[i]) && partsColl[i].mother() &&
        std::abs(partsColl[i].mother()->pdgId()) == TopDecayID::WID) {
      cand = &partsColl[i];
    }
  }
  return cand;
}

const reco::GenParticle* StGenEvent::singleW() const {
  const reco::GenParticle* cand = nullptr;
  if (singleLepton()) {
    const reco::GenParticleCollection& partsColl = *parts_;
    const reco::GenParticle& singleLep = *singleLepton();
    for (unsigned int i = 0; i < partsColl.size(); ++i) {
      if (std::abs(partsColl[i].pdgId()) == TopDecayID::WID &&
          reco::flavour(singleLep) == -reco::flavour(partsColl[i])) {
        // PDG Id:13=mu- 24=W+ (+24)->(-13) (-24)->(+13) opposite sign
        cand = &partsColl[i];
      }
    }
  }
  return cand;
}

const reco::GenParticle* StGenEvent::singleTop() const {
  const reco::GenParticle* cand = nullptr;
  if (singleLepton()) {
    const reco::GenParticleCollection& partsColl = *parts_;
    const reco::GenParticle& singleLep = *singleLepton();
    for (unsigned int i = 0; i < partsColl.size(); ++i) {
      if (std::abs(partsColl[i].pdgId()) == TopDecayID::tID &&
          reco::flavour(singleLep) != reco::flavour(partsColl[i])) {
        cand = &partsColl[i];
      }
    }
  }
  return cand;
}