#include <TtSemiLepHypHitFit.h>
Public Member Functions | |
TtSemiLepHypHitFit (const edm::ParameterSet &) | |
~TtSemiLepHypHitFit () | |
Private Member Functions | |
virtual 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) |
build event hypothesis from the reco objects of a semi-leptonic event | |
virtual void | buildKey () |
build the event hypothesis key | |
Private Attributes | |
edm::InputTag | leptons_ |
edm::InputTag | neutrinos_ |
edm::InputTag | partonsHadB_ |
edm::InputTag | partonsHadP_ |
edm::InputTag | partonsHadQ_ |
edm::InputTag | partonsLepB_ |
edm::InputTag | status_ |
Definition at line 6 of file TtSemiLepHypHitFit.h.
TtSemiLepHypHitFit::TtSemiLepHypHitFit | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 5 of file TtSemiLepHypHitFit.cc.
: TtSemiLepHypothesis( cfg ), status_ (cfg.getParameter<edm::InputTag>("status" )), partonsHadP_(cfg.getParameter<edm::InputTag>("partonsHadP")), partonsHadQ_(cfg.getParameter<edm::InputTag>("partonsHadQ")), partonsHadB_(cfg.getParameter<edm::InputTag>("partonsHadB")), partonsLepB_(cfg.getParameter<edm::InputTag>("partonsLepB")), leptons_ (cfg.getParameter<edm::InputTag>("leptons" )), neutrinos_ (cfg.getParameter<edm::InputTag>("neutrinos" )) { }
TtSemiLepHypHitFit::~TtSemiLepHypHitFit | ( | ) |
Definition at line 17 of file TtSemiLepHypHitFit.cc.
{ }
void TtSemiLepHypHitFit::buildHypo | ( | edm::Event & | evt, |
const edm::Handle< edm::View< reco::RecoCandidate > > & | leps, | ||
const edm::Handle< std::vector< pat::MET > > & | mets, | ||
const edm::Handle< std::vector< pat::Jet > > & | jets, | ||
std::vector< int > & | match, | ||
const unsigned int | iComb | ||
) | [private, virtual] |
build event hypothesis from the reco objects of a semi-leptonic event
Implements TtSemiLepHypothesis.
Definition at line 20 of file TtSemiLepHypHitFit.cc.
References edm::Event::getByLabel(), TtSemiLepHypothesis::hadronicB_, TtSemiLepHypothesis::lepton_, TtSemiLepHypothesis::leptonicB_, EgammaValidation_Wenu_cff::leptons, leptons_, TtSemiLepHypothesis::lightQ_, TtSemiLepHypothesis::lightQBar_, TtSemiLepHypothesis::neutrino_, neutrinos_, partonsHadB_, partonsHadP_, partonsHadQ_, partonsLepB_, TtSemiLepHypothesis::setCandidate(), ntuplemaker::status, and status_.
{ edm::Handle<std::vector<int> > status; evt.getByLabel(status_, status); if( (*status)[iComb] != 0 ){ // create empty hypothesis if kinematic fit did not converge return; } edm::Handle<std::vector<pat::Particle> > partonsHadP; edm::Handle<std::vector<pat::Particle> > partonsHadQ; edm::Handle<std::vector<pat::Particle> > partonsHadB; edm::Handle<std::vector<pat::Particle> > partonsLepB; edm::Handle<std::vector<pat::Particle> > leptons; edm::Handle<std::vector<pat::Particle> > neutrinos; evt.getByLabel(partonsHadP_, partonsHadP); evt.getByLabel(partonsHadQ_, partonsHadQ); evt.getByLabel(partonsHadB_, partonsHadB); evt.getByLabel(partonsLepB_, partonsLepB); evt.getByLabel(leptons_ , leptons ); evt.getByLabel(neutrinos_ , neutrinos ); // ----------------------------------------------------- // add jets // ----------------------------------------------------- if( !( partonsHadP->empty() || partonsHadQ->empty() || partonsHadB->empty() || partonsLepB->empty() ) ) { setCandidate(partonsHadP, iComb, lightQ_ ); setCandidate(partonsHadQ, iComb, lightQBar_); setCandidate(partonsHadB, iComb, hadronicB_); setCandidate(partonsLepB, iComb, leptonicB_); } // ----------------------------------------------------- // add lepton // ----------------------------------------------------- if( !leptons->empty() ){ setCandidate(leptons, iComb, lepton_); } match.push_back( 0 ); // ----------------------------------------------------- // add neutrino // ----------------------------------------------------- if( !neutrinos->empty() ){ setCandidate(neutrinos, iComb, neutrino_); } }
virtual void TtSemiLepHypHitFit::buildKey | ( | ) | [inline, private, virtual] |
build the event hypothesis key
Implements TtSemiLepHypothesis.
Definition at line 16 of file TtSemiLepHypHitFit.h.
References TtSemiLepHypothesis::key_, and TtEvent::kHitFit.
{ key_= TtSemiLeptonicEvent::kHitFit; };
edm::InputTag TtSemiLepHypHitFit::leptons_ [private] |
Definition at line 29 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
edm::InputTag TtSemiLepHypHitFit::neutrinos_ [private] |
Definition at line 30 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
Definition at line 27 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
Definition at line 25 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
Definition at line 26 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
Definition at line 28 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().
edm::InputTag TtSemiLepHypHitFit::status_ [private] |
Definition at line 24 of file TtSemiLepHypHitFit.h.
Referenced by buildHypo().