CMS 3D CMS Logo

TtSemiLepHypGenMatch.cc
Go to the documentation of this file.
5 
7  : TtSemiLepHypothesis(cfg), genEvtToken_(consumes<TtGenEvent>(edm::InputTag("genEvt"))) {}
8 
10 
13  const edm::Handle<std::vector<pat::MET> >& mets,
14  const edm::Handle<std::vector<pat::Jet> >& jets,
15  std::vector<int>& match,
16  const unsigned int iComb) {
17  // -----------------------------------------------------
18  // get genEvent (to distinguish between uds and c quarks
19  // and for the lepton matching)
20  // -----------------------------------------------------
22  evt.getByToken(genEvtToken_, genEvt);
23 
24  // -----------------------------------------------------
25  // add jets
26  // -----------------------------------------------------
27  for (unsigned idx = 0; idx < match.size(); ++idx) {
28  if (isValid(match[idx], jets)) {
29  switch (idx) {
31  if (std::abs(genEvt->hadronicDecayQuark()->pdgId()) == 4)
32  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("cQuark"));
33  else
34  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("udsQuark"));
35  break;
37  if (std::abs(genEvt->hadronicDecayQuarkBar()->pdgId()) == 4)
38  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("cQuark"));
39  else
40  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("udsQuark"));
41  break;
43  setCandidate(jets, match[idx], hadronicB_, jetCorrectionLevel("bQuark"));
44  break;
46  setCandidate(jets, match[idx], leptonicB_, jetCorrectionLevel("bQuark"));
47  break;
48  }
49  }
50  }
51 
52  // -----------------------------------------------------
53  // add lepton
54  // -----------------------------------------------------
55  int iLepton = findMatchingLepton(genEvt, leps);
56  if (iLepton < 0)
57  return;
58  setCandidate(leps, iLepton, lepton_);
59  match.push_back(iLepton);
60 
61  // -----------------------------------------------------
62  // add neutrino
63  // -----------------------------------------------------
64  if (mets->empty())
65  return;
66  if (neutrinoSolutionType_ == -1)
68  else
70 }
71 
75  int genIdx = -1;
76 
77  // jump out with -1 when the collection is empty
78  if (leps->empty())
79  return genIdx;
80 
81  if (genEvt->isTtBar() && genEvt->isSemiLeptonic(leptonType(&(leps->front()))) && genEvt->singleLepton()) {
82  double minDR = -1;
83  for (unsigned i = 0; i < leps->size(); ++i) {
84  double dR =
85  deltaR(genEvt->singleLepton()->eta(), genEvt->singleLepton()->phi(), (*leps)[i].eta(), (*leps)[i].phi());
86  if (minDR < 0 || dR < minDR) {
87  minDR = dR;
88  genIdx = i;
89  }
90  }
91  }
92  return genIdx;
93 }
int pdgId() const final
PDG identifier.
bool isSemiLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as semi-laptonic
Definition: TtGenEvent.h:38
double eta() const final
momentum pseudorapidity
WDecay::LeptonType leptonType(const reco::RecoCandidate *cand)
determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kN...
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
reco::ShallowClonePtrCandidate * lepton_
const reco::GenParticle * hadronicDecayQuarkBar() const
get light anti-quark of hadronic decay branch
Definition: TtGenEvent.h:78
reco::ShallowClonePtrCandidate * lightQBar_
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
reco::ShallowClonePtrCandidate * neutrino_
void buildHypo(edm::Event &, const edm::Handle< edm::View< reco::RecoCandidate > > &, const edm::Handle< std::vector< pat::MET > > &, const edm::Handle< std::vector< pat::Jet > > &, std::vector< int > &, const unsigned int iComb) override
build event hypothesis from the reco objects of a semi-leptonic event
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const reco::GenParticle * hadronicDecayQuark(bool invertFlavor=false) const
get light quark of hadronic decay branch
Definition: TtGenEvent.cc:153
int findMatchingLepton(const edm::Handle< TtGenEvent > &genEvt, const edm::Handle< edm::View< reco::RecoCandidate > > &)
find index of the candidate nearest to the singleLepton of the generator event in the collection; ret...
reco::ShallowClonePtrCandidate * hadronicB_
TtSemiLepHypGenMatch(const edm::ParameterSet &)
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const int &idx, const int &type)
set neutrino, using mW = 80.4 to calculate the neutrino pz
HLT enums.
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
int neutrinoSolutionType_
algorithm used to calculate neutrino solutions (see cfi for further details)
edm::EDGetTokenT< TtGenEvent > genEvtToken_
double phi() const final
momentum azimuthal angle
reco::ShallowClonePtrCandidate * lightQ_
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:28
reco::ShallowClonePtrCandidate * leptonicB_
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
const reco::GenParticle * singleLepton(bool excludeTauLeptons=false) const
return single lepton if available; 0 else
Definition: TtGenEvent.cc:97