CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullHadHypGenMatch.cc
Go to the documentation of this file.
3 
5  : TtFullHadHypothesis(cfg), genEvtToken_(consumes<TtGenEvent>(edm::InputTag("genEvt"))) {}
6 
8 
10  const edm::Handle<std::vector<pat::Jet> >& jets,
11  std::vector<int>& match,
12  const unsigned int iComb) {
13  // -----------------------------------------------------
14  // get genEvent (to distinguish between uds and c quarks)
15  // -----------------------------------------------------
17  evt.getByToken(genEvtToken_, genEvt);
18 
19  // -----------------------------------------------------
20  // add jets
21  // -----------------------------------------------------
22  for (unsigned idx = 0; idx < match.size(); ++idx) {
23  if (isValid(match[idx], jets)) {
24  switch (idx) {
26  if (std::abs(genEvt->daughterQuarkOfWPlus()->pdgId()) == 4)
27  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("cQuark"));
28  else
29  setCandidate(jets, match[idx], lightQ_, jetCorrectionLevel("udsQuark"));
30  break;
32  if (std::abs(genEvt->daughterQuarkBarOfWPlus()->pdgId()) == 4)
33  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("cQuark"));
34  else
35  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("udsQuark"));
36  break;
38  setCandidate(jets, match[idx], b_, jetCorrectionLevel("bQuark"));
39  break;
41  if (std::abs(genEvt->daughterQuarkOfWMinus()->pdgId()) == 4)
42  setCandidate(jets, match[idx], lightP_, jetCorrectionLevel("cQuark"));
43  else
44  setCandidate(jets, match[idx], lightP_, jetCorrectionLevel("udsQuark"));
45  break;
47  if (std::abs(genEvt->daughterQuarkBarOfWMinus()->pdgId()) == 4)
48  setCandidate(jets, match[idx], lightPBar_, jetCorrectionLevel("cQuark"));
49  else
50  setCandidate(jets, match[idx], lightPBar_, jetCorrectionLevel("udsQuark"));
51  break;
53  setCandidate(jets, match[idx], bBar_, jetCorrectionLevel("bQuark"));
54  break;
55  }
56  }
57  }
58 }
tuple cfg
Definition: looper.py:296
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
reco::ShallowClonePtrCandidate * lightPBar_
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
vector< PseudoJet > jets
reco::ShallowClonePtrCandidate * lightQBar_
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::ShallowClonePtrCandidate * lightP_
void buildHypo(edm::Event &evt, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &match, const unsigned int iComb) override
build event hypothesis from the reco objects of a semi-leptonic event
reco::ShallowClonePtrCandidate * bBar_
reco::ShallowClonePtrCandidate * b_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
edm::EDGetTokenT< TtGenEvent > genEvtToken_
genEvtToken_(mayConsume< TtGenEvent >(genEvt_))
reco::ShallowClonePtrCandidate * lightQ_
TtFullHadHypGenMatch(const edm::ParameterSet &cfg)