CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullHadHypGenMatch.cc
Go to the documentation of this file.
4 
6  TtFullHadHypothesis( cfg )
7 {
8 }
9 
11 
12 void
14  const edm::Handle<std::vector<pat::Jet> >& jets,
15  std::vector<int>& match,
16  const unsigned int iComb)
17 {
18  // -----------------------------------------------------
19  // get genEvent (to distinguish between uds and c quarks
20  // and for the lepton matching)
21  // -----------------------------------------------------
23  evt.getByLabel("genEvt", genEvt);
24 
25  // -----------------------------------------------------
26  // add jets
27  // -----------------------------------------------------
28  for(unsigned idx=0; idx<match.size(); ++idx){
29  if( isValid(match[idx], jets) ){
30  switch(idx){
32  if( std::abs(genEvt->daughterQuarkOfWPlus()->pdgId())==4 )
33  setCandidate(jets, match[idx], lightQ_ , jetCorrectionLevel("cQuark"));
34  else
35  setCandidate(jets, match[idx], lightQ_ , jetCorrectionLevel("udsQuark"));
36  break;
38  if( std::abs(genEvt->daughterQuarkBarOfWPlus()->pdgId())==4 )
39  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("cQuark"));
40  else
41  setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("udsQuark"));
42  break;
44  setCandidate(jets, match[idx], b_ , jetCorrectionLevel("bQuark")); break;
46  if( std::abs(genEvt->daughterQuarkOfWMinus()->pdgId())==4 )
47  setCandidate(jets, match[idx], lightP_ , jetCorrectionLevel("cQuark"));
48  else
49  setCandidate(jets, match[idx], lightP_ , jetCorrectionLevel("udsQuark"));
50  break;
52  if( std::abs(genEvt->daughterQuarkBarOfWMinus()->pdgId())==4 )
53  setCandidate(jets, match[idx], lightPBar_, jetCorrectionLevel("cQuark"));
54  else
55  setCandidate(jets, match[idx], lightPBar_, jetCorrectionLevel("udsQuark"));
56  break;
58  setCandidate(jets, match[idx], bBar_ , jetCorrectionLevel("bQuark")); break;
59  }
60  }
61  }
62 }
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
#define abs(x)
Definition: mlp_lapack.h:159
reco::ShallowClonePtrCandidate * lightPBar_
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
reco::ShallowClonePtrCandidate * lightP_
virtual void buildHypo(edm::Event &evt, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &match, const unsigned int iComb)
build event hypothesis from the reco objects of a semi-leptonic event
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
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:6
reco::ShallowClonePtrCandidate * lightQ_
TtFullHadHypGenMatch(const edm::ParameterSet &cfg)