CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullLepHypKinSolution.cc
Go to the documentation of this file.
3 
4 
6  TtFullLepHypothesis( cfg ),
7  nusToken_ (consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("Neutrinos" ))),
8  nuBarsToken_ (consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("NeutrinoBars" ))),
9  solWeightToken_ (consumes<std::vector<double> >(cfg.getParameter<edm::InputTag>("solutionWeight" )))
10 {
11 }
12 
14 
15 void
17  const edm::Handle<std::vector<pat::Electron > >& elecs,
18  const edm::Handle<std::vector<pat::Muon> >& mus,
19  const edm::Handle<std::vector<pat::Jet> >& jets,
20  const edm::Handle<std::vector<pat::MET> >& mets,
21  std::vector<int>& match,
22  const unsigned int iComb)
23 {
25 // edm::Handle<std::vector<std::vector<int> > > idcsVec;
28 
29  evt.getByToken(solWeightToken_, solWeight);
30 // evt.getByToken(particleIdcsToken_, idcsVec );
31  evt.getByToken(nusToken_, nus );
32  evt.getByToken(nuBarsToken_, nuBars );
33 
34  if( (*solWeight)[iComb]<0 ){
35  // create empty hypothesis if no solution exists
36  return;
37  }
38 
39  // -----------------------------------------------------
40  // add jets
41  // -----------------------------------------------------
42  if( !jets->empty() ){
45  }
46  // -----------------------------------------------------
47  // add leptons
48  // -----------------------------------------------------
49  if( !elecs->empty() && match[2]>=0)
50  setCandidate(elecs, match[2], leptonBar_);
51 
52  if( !elecs->empty() && match[3]>=0)
53  setCandidate(elecs, match[3], lepton_);
54 
55  if( !mus->empty() && match[4]>=0 && match[2]<0)
56  setCandidate(mus, match[4], leptonBar_);
57 
58  // this 'else' happens if you have a wrong charge electron-muon-
59  // solution so the indices are (b-idx, bbar-idx, 0, -1, 0, -1)
60  // so the mu^+ is stored as l^-
61  else if( !mus->empty() && match[4]>=0)
62  setCandidate(mus, match[4], lepton_);
63 
64  if( !mus->empty() && match[5]>=0 && match[3]<0)
65  setCandidate(mus, match[5], lepton_);
66 
67  // this 'else' happens if you have a wrong charge electron-muon-
68  // solution so the indices are (b-idx, bbar-idx, -1, 0, -1, 0)
69  // so the mu^- is stored as l^+
70  else if( !mus->empty() && match[5]>=0)
71  setCandidate(mus, match[5], leptonBar_);
72 
73  // -----------------------------------------------------
74  // add neutrinos
75  // -----------------------------------------------------
76  if( !nus->empty() )
77  setCandidate(nus, iComb, neutrino_);
78 
79  if( !nuBars->empty() )
80  setCandidate(nuBars, iComb, neutrinoBar_);
81 
82 }
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * leptonBar_
tuple cfg
Definition: looper.py:293
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:464
edm::EDGetTokenT< std::vector< reco::LeafCandidate > > nuBarsToken_
virtual 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)
build event hypothesis from the reco objects of a full-leptonic event
vector< PseudoJet > jets
reco::ShallowClonePtrCandidate * lepton_
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_