CMS 3D CMS Logo

TtSemiLepHypHitFit.cc
Go to the documentation of this file.
2 
4  : TtSemiLepHypothesis(cfg),
5  statusToken_(consumes<std::vector<int> >(cfg.getParameter<edm::InputTag>("status"))),
6  partonsHadPToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadP"))),
7  partonsHadQToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadQ"))),
8  partonsHadBToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadB"))),
9  partonsLepBToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsLepB"))),
10  leptonsToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("leptons"))),
11  neutrinosToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("neutrinos"))) {}
12 
14 
17  const edm::Handle<std::vector<pat::MET> >& mets,
18  const edm::Handle<std::vector<pat::Jet> >& jets,
19  std::vector<int>& match,
20  const unsigned int iComb) {
22  evt.getByToken(statusToken_, status);
23  if ((*status)[iComb] != 0) {
24  // create empty hypothesis if kinematic fit did not converge
25  return;
26  }
27 
34 
35  evt.getByToken(partonsHadPToken_, partonsHadP);
36  evt.getByToken(partonsHadQToken_, partonsHadQ);
37  evt.getByToken(partonsHadBToken_, partonsHadB);
38  evt.getByToken(partonsLepBToken_, partonsLepB);
39  evt.getByToken(leptonsToken_, leptons);
40  evt.getByToken(neutrinosToken_, neutrinos);
41 
42  // -----------------------------------------------------
43  // add jets
44  // -----------------------------------------------------
45  if (!(partonsHadP->empty() || partonsHadQ->empty() || partonsHadB->empty() || partonsLepB->empty())) {
46  setCandidate(partonsHadP, iComb, lightQ_);
47  setCandidate(partonsHadQ, iComb, lightQBar_);
48  setCandidate(partonsHadB, iComb, hadronicB_);
49  setCandidate(partonsLepB, iComb, leptonicB_);
50  }
51 
52  // -----------------------------------------------------
53  // add lepton
54  // -----------------------------------------------------
55  if (!leptons->empty()) {
56  setCandidate(leptons, iComb, lepton_);
57  }
58  match.push_back(0);
59 
60  // -----------------------------------------------------
61  // add neutrino
62  // -----------------------------------------------------
63  if (!neutrinos->empty()) {
64  setCandidate(neutrinos, iComb, neutrino_);
65  }
66 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadPToken_
~TtSemiLepHypHitFit() override
edm::EDGetTokenT< std::vector< int > > statusToken_
edm::EDGetTokenT< std::vector< pat::Particle > > leptonsToken_
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
TtSemiLepHypHitFit(const edm::ParameterSet &)
Definition: HeavyIon.h:7
reco::ShallowClonePtrCandidate * neutrino_
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadQToken_
reco::ShallowClonePtrCandidate * hadronicB_
edm::EDGetTokenT< std::vector< pat::Particle > > partonsLepBToken_
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadBToken_
HLT enums.
edm::EDGetTokenT< std::vector< pat::Particle > > neutrinosToken_
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
reco::ShallowClonePtrCandidate * lightQ_
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
reco::ShallowClonePtrCandidate * leptonicB_