CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullLepHypKinSolution.cc
Go to the documentation of this file.
3 
5  : TtFullLepHypothesis(cfg),
6  nusToken_(consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("Neutrinos"))),
7  nuBarsToken_(consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("NeutrinoBars"))),
8  solWeightToken_(consumes<std::vector<double> >(cfg.getParameter<edm::InputTag>("solutionWeight"))) {}
9 
11 
13  const edm::Handle<std::vector<pat::Electron> >& elecs,
14  const edm::Handle<std::vector<pat::Muon> >& mus,
15  const edm::Handle<std::vector<pat::Jet> >& jets,
16  const edm::Handle<std::vector<pat::MET> >& mets,
17  std::vector<int>& match,
18  const unsigned int iComb) {
20  // edm::Handle<std::vector<std::vector<int> > > idcsVec;
23 
24  evt.getByToken(solWeightToken_, solWeight);
25  // evt.getByToken(particleIdcsToken_, idcsVec );
26  evt.getByToken(nusToken_, nus);
27  evt.getByToken(nuBarsToken_, nuBars);
28 
29  if ((*solWeight)[iComb] < 0) {
30  // create empty hypothesis if no solution exists
31  return;
32  }
33 
34  // -----------------------------------------------------
35  // add jets
36  // -----------------------------------------------------
37  if (!jets->empty()) {
40  }
41  // -----------------------------------------------------
42  // add leptons
43  // -----------------------------------------------------
44  if (!elecs->empty() && match[2] >= 0)
45  setCandidate(elecs, match[2], leptonBar_);
46 
47  if (!elecs->empty() && match[3] >= 0)
48  setCandidate(elecs, match[3], lepton_);
49 
50  if (!mus->empty() && match[4] >= 0 && match[2] < 0)
51  setCandidate(mus, match[4], leptonBar_);
52 
53  // this 'else' happens if you have a wrong charge electron-muon-
54  // solution so the indices are (b-idx, bbar-idx, 0, -1, 0, -1)
55  // so the mu^+ is stored as l^-
56  else if (!mus->empty() && match[4] >= 0)
57  setCandidate(mus, match[4], lepton_);
58 
59  if (!mus->empty() && match[5] >= 0 && match[3] < 0)
60  setCandidate(mus, match[5], lepton_);
61 
62  // this 'else' happens if you have a wrong charge electron-muon-
63  // solution so the indices are (b-idx, bbar-idx, -1, 0, -1, 0)
64  // so the mu^- is stored as l^+
65  else if (!mus->empty() && match[5] >= 0)
66  setCandidate(mus, match[5], leptonBar_);
67 
68  // -----------------------------------------------------
69  // add neutrinos
70  // -----------------------------------------------------
71  if (!nus->empty())
72  setCandidate(nus, iComb, neutrino_);
73 
74  if (!nuBars->empty())
75  setCandidate(nuBars, iComb, neutrinoBar_);
76 }
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * leptonBar_
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
edm::EDGetTokenT< std::vector< reco::LeafCandidate > > nuBarsToken_
vector< PseudoJet > jets
reco::ShallowClonePtrCandidate * lepton_
void buildHypo(edm::Event &evt, const edm::Handle< std::vector< pat::Electron > > &elecs, const edm::Handle< std::vector< pat::Muon > > &mus, const edm::Handle< std::vector< pat::Jet > > &jets, const edm::Handle< std::vector< pat::MET > > &mets, std::vector< int > &match, const unsigned int iComb) override
build event hypothesis from the reco objects of a full-leptonic event
reco::ShallowClonePtrCandidate * bBar_
std::string jetCorrectionLevel_
edm::EDGetTokenT< std::vector< reco::LeafCandidate > > nusToken_
edm::EDGetTokenT< std::vector< double > > solWeightToken_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
TtFullLepHypKinSolution(const edm::ParameterSet &)
reco::ShallowClonePtrCandidate * b_
reco::ShallowClonePtrCandidate * neutrinoBar_