CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtSemiLepHypHitFit.cc
Go to the documentation of this file.
2 
4 
6 public:
7  explicit TtSemiLepHypHitFit(const edm::ParameterSet&);
8 
9 private:
11  void buildKey() override { key_ = TtSemiLeptonicEvent::kHitFit; };
13  void buildHypo(edm::Event&,
15  const edm::Handle<std::vector<pat::MET> >&,
16  const edm::Handle<std::vector<pat::Jet> >&,
17  std::vector<int>&,
18  const unsigned int iComb) override;
19 
27 };
28 
30  : TtSemiLepHypothesis(cfg),
31  statusToken_(consumes<std::vector<int> >(cfg.getParameter<edm::InputTag>("status"))),
32  partonsHadPToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadP"))),
33  partonsHadQToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadQ"))),
34  partonsHadBToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsHadB"))),
35  partonsLepBToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("partonsLepB"))),
36  leptonsToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("leptons"))),
37  neutrinosToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("neutrinos"))) {}
38 
41  const edm::Handle<std::vector<pat::MET> >& mets,
42  const edm::Handle<std::vector<pat::Jet> >& jets,
43  std::vector<int>& match,
44  const unsigned int iComb) {
46  evt.getByToken(statusToken_, status);
47  if ((*status)[iComb] != 0) {
48  // create empty hypothesis if kinematic fit did not converge
49  return;
50  }
51 
58 
59  evt.getByToken(partonsHadPToken_, partonsHadP);
60  evt.getByToken(partonsHadQToken_, partonsHadQ);
61  evt.getByToken(partonsHadBToken_, partonsHadB);
62  evt.getByToken(partonsLepBToken_, partonsLepB);
63  evt.getByToken(leptonsToken_, leptons);
64  evt.getByToken(neutrinosToken_, neutrinos);
65 
66  // -----------------------------------------------------
67  // add jets
68  // -----------------------------------------------------
69  if (!(partonsHadP->empty() || partonsHadQ->empty() || partonsHadB->empty() || partonsLepB->empty())) {
70  lightQ_ = makeCandidate(partonsHadP, iComb);
71  lightQBar_ = makeCandidate(partonsHadQ, iComb);
72  hadronicB_ = makeCandidate(partonsHadB, iComb);
73  leptonicB_ = makeCandidate(partonsLepB, iComb);
74  }
75 
76  // -----------------------------------------------------
77  // add lepton
78  // -----------------------------------------------------
79  if (!leptons->empty()) {
80  lepton_ = makeCandidate(leptons, iComb);
81  }
82  match.push_back(0);
83 
84  // -----------------------------------------------------
85  // add neutrino
86  // -----------------------------------------------------
87  if (!neutrinos->empty()) {
88  neutrino_ = makeCandidate(neutrinos, iComb);
89  }
90 }
91 
tuple cfg
Definition: looper.py:296
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadPToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< int > > statusToken_
edm::EDGetTokenT< std::vector< pat::Particle > > leptonsToken_
list status
Definition: mps_update.py:107
TtSemiLepHypHitFit(const edm::ParameterSet &)
void buildKey() override
build the event hypothesis key
std::unique_ptr< reco::ShallowClonePtrCandidate > leptonicB_
std::unique_ptr< reco::ShallowClonePtrCandidate > lepton_
vector< PseudoJet > jets
std::unique_ptr< reco::ShallowClonePtrCandidate > makeCandidate(const edm::Handle< C > &handle, const int &idx)
use one object in a collection to set a ShallowClonePtrCandidate
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadQToken_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQBar_
edm::EDGetTokenT< std::vector< pat::Particle > > partonsLepBToken_
edm::EDGetTokenT< std::vector< pat::Particle > > partonsHadBToken_
std::unique_ptr< reco::ShallowClonePtrCandidate > neutrino_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQ_
edm::EDGetTokenT< std::vector< pat::Particle > > neutrinosToken_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
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
std::unique_ptr< reco::ShallowClonePtrCandidate > hadronicB_
int key_
hypothesis key (to be set by the buildKey function)