CMS 3D CMS Logo

TtSemiLepHypothesis Class Reference

#include <TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepHypothesis.h>

Inheritance diagram for TtSemiLepHypothesis:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper TtSemiLepHypGenMatch TtSemiLepHypGeom TtSemiLepHypKinFit TtSemiLepHypMaxSumPtWMass TtSemiLepHypMVADisc TtSemiLepHypWMassMaxSumPt

List of all members.

Public Member Functions

 TtSemiLepHypothesis (const edm::ParameterSet &)
 ~TtSemiLepHypothesis ()

Protected Member Functions

virtual void buildHypo (edm::Event &event, const edm::Handle< edm::View< reco::RecoCandidate > > &lepton, const edm::Handle< std::vector< pat::MET > > &neutrino, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation, const unsigned int iComb)=0
 build event hypothesis from the reco objects of a semi-leptonic event
virtual void buildKey ()=0
 build the event hypothesis key
reco::CompositeCandidate hypo ()
 return event hypothesis
bool isValid (const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
 check if index is in valid range of selected jets
int key () const
 return key
virtual void produce (edm::Event &, const edm::EventSetup &)
 produce the event hypothesis as CompositeCandidate and Key
void resetCandidates ()
 reset candidate pointers before hypo build process
template<typename C>
void setCandidate (const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
 use one object in a collection to set a ShallowClonePtrCandidate

Protected Attributes

bool getMatch_
reco::ShallowClonePtrCandidatehadronicB_
edm::InputTag jets_
int key_
edm::InputTag leps_
reco::ShallowClonePtrCandidatelepton_
reco::ShallowClonePtrCandidateleptonicB_
reco::ShallowClonePtrCandidatelightQ_
reco::ShallowClonePtrCandidatelightQBar_
edm::InputTag match_
edm::InputTag mets_
reco::ShallowClonePtrCandidateneutrino_


Detailed Description

Definition at line 19 of file TtSemiLepHypothesis.h.


Constructor & Destructor Documentation

TtSemiLepHypothesis::TtSemiLepHypothesis ( const edm::ParameterSet cfg  )  [explicit]

Definition at line 5 of file TtSemiLepHypothesis.cc.

References edm::ParameterSet::exists(), getMatch_, edm::ParameterSet::getParameter(), and match_.

00005                                                                   :
00006   jets_(cfg.getParameter<edm::InputTag>("jets")),
00007   leps_(cfg.getParameter<edm::InputTag>("leps")),
00008   mets_(cfg.getParameter<edm::InputTag>("mets")),
00009   lightQ_(0), lightQBar_(0), hadronicB_(0), 
00010   leptonicB_(0), neutrino_(0), lepton_(0)
00011 {
00012   getMatch_ = false;
00013   if( cfg.exists("match") ) {
00014     getMatch_ = true;
00015     match_ = cfg.getParameter<edm::InputTag>("match");
00016   }
00017 
00018   produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
00019   produces<int>("Key");
00020 }

TtSemiLepHypothesis::~TtSemiLepHypothesis (  ) 

Definition at line 22 of file TtSemiLepHypothesis.cc.

References hadronicB_, lepton_, leptonicB_, lightQ_, lightQBar_, and neutrino_.

00023 {
00024   if( lightQ_    ) delete lightQ_;
00025   if( lightQBar_ ) delete lightQBar_;
00026   if( hadronicB_ ) delete hadronicB_;
00027   if( leptonicB_ ) delete leptonicB_;
00028   if( neutrino_  ) delete neutrino_;
00029   if( lepton_    ) delete lepton_;
00030 }


Member Function Documentation

virtual void TtSemiLepHypothesis::buildHypo ( edm::Event event,
const edm::Handle< edm::View< reco::RecoCandidate > > &  lepton,
const edm::Handle< std::vector< pat::MET > > &  neutrino,
const edm::Handle< std::vector< pat::Jet > > &  jets,
std::vector< int > &  jetPartonAssociation,
const unsigned int  iComb 
) [protected, pure virtual]

build event hypothesis from the reco objects of a semi-leptonic event

Implemented in TtSemiLepHypGenMatch, TtSemiLepHypGeom, TtSemiLepHypKinFit, TtSemiLepHypMaxSumPtWMass, TtSemiLepHypMVADisc, and TtSemiLepHypWMassMaxSumPt.

Referenced by produce().

virtual void TtSemiLepHypothesis::buildKey (  )  [protected, pure virtual]

build the event hypothesis key

Implemented in TtSemiLepHypGenMatch, TtSemiLepHypGeom, TtSemiLepHypKinFit, TtSemiLepHypMaxSumPtWMass, TtSemiLepHypMVADisc, and TtSemiLepHypWMassMaxSumPt.

Referenced by produce().

reco::CompositeCandidate TtSemiLepHypothesis::hypo (  )  [protected]

return event hypothesis

Definition at line 97 of file TtSemiLepHypothesis.cc.

References reco::CompositeCandidate::addDaughter(), TtSemiLepDaughter::HadB, TtSemiLepDaughter::HadP, TtSemiLepDaughter::HadQ, hadronicB_, TtSemiLepDaughter::HadTop, TtSemiLepDaughter::HadW, TtSemiLepDaughter::Lep, TtSemiLepDaughter::LepB, lepton_, leptonicB_, TtSemiLepDaughter::LepTop, TtSemiLepDaughter::LepW, lightQ_, lightQBar_, neutrino_, TtSemiLepDaughter::Nu, and AddFourMomenta::set().

Referenced by produce().

00098 {
00099   // check for sanity of the hypothesis
00100   if( !lightQ_ || !lightQBar_ || !hadronicB_ || 
00101       !leptonicB_ || !neutrino_ || !lepton_ )
00102     return reco::CompositeCandidate();
00103   
00104   // setup transient references
00105   reco::CompositeCandidate hyp, hadTop, hadW, lepTop, lepW;
00106 
00107   AddFourMomenta addFourMomenta;  
00108   // build up the top branch that decays leptonically
00109   lepW  .addDaughter(*lepton_,   TtSemiLepDaughter::Lep    );
00110   lepW  .addDaughter(*neutrino_, TtSemiLepDaughter::Nu     );
00111   addFourMomenta.set( lepW );
00112   lepTop.addDaughter( lepW,      TtSemiLepDaughter::LepW   );
00113   lepTop.addDaughter(*leptonicB_,TtSemiLepDaughter::LepB   );
00114   addFourMomenta.set( lepTop );
00115   
00116   // build up the top branch that decays hadronically
00117   hadW  .addDaughter(*lightQ_,   TtSemiLepDaughter::HadP   );
00118   hadW  .addDaughter(*lightQBar_,TtSemiLepDaughter::HadQ   );
00119   addFourMomenta.set( hadW );
00120   hadTop.addDaughter( hadW,      TtSemiLepDaughter::HadW   );
00121   hadTop.addDaughter(*hadronicB_,TtSemiLepDaughter::HadB   );
00122   addFourMomenta.set( hadTop );
00123 
00124   // build ttbar hypotheses
00125   hyp.addDaughter( lepTop,       TtSemiLepDaughter::LepTop );
00126   hyp.addDaughter( hadTop,       TtSemiLepDaughter::HadTop );
00127   addFourMomenta.set( hyp );
00128 
00129   return hyp;
00130 }

bool TtSemiLepHypothesis::isValid ( const int idx,
const edm::Handle< std::vector< pat::Jet > > &  jets 
) [inline, protected]

check if index is in valid range of selected jets

Definition at line 40 of file TtSemiLepHypothesis.h.

References pfTauBenchmarkGeneric_cfi::jets.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), and TtSemiLepHypMVADisc::buildHypo().

00040 { return (0<=idx && idx<(int)jets->size()); };

int TtSemiLepHypothesis::key (  )  const [inline, protected]

return key

Definition at line 36 of file TtSemiLepHypothesis.h.

References key_.

Referenced by produce().

00036 { return key_; };

void TtSemiLepHypothesis::produce ( edm::Event evt,
const edm::EventSetup setup 
) [protected, virtual]

produce the event hypothesis as CompositeCandidate and Key

Implements edm::EDProducer.

Definition at line 33 of file TtSemiLepHypothesis.cc.

References buildHypo(), buildKey(), edm::Event::getByLabel(), getMatch_, hypo(), i, pfTauBenchmarkGeneric_cfi::jets, jets_, key(), leps_, edm::match(), match_, mets_, edm::Event::put(), and resetCandidates().

00034 {
00035   edm::Handle<std::vector<pat::Jet> > jets;
00036   evt.getByLabel(jets_, jets);
00037   
00038   edm::Handle<edm::View<reco::RecoCandidate> > leps;
00039   evt.getByLabel(leps_, leps);
00040 
00041   edm::Handle<std::vector<pat::MET> > mets;
00042   evt.getByLabel(mets_, mets);
00043 
00044   std::vector<std::vector<int> > matchVec;
00045   if( getMatch_ ) {
00046     edm::Handle<std::vector<std::vector<int> > > matchHandle;
00047     evt.getByLabel(match_, matchHandle);
00048     matchVec = *matchHandle;
00049   }
00050   else {
00051     std::vector<int> dummyMatch;
00052     for(unsigned int i = 0; i < 4; ++i) 
00053       dummyMatch.push_back( -1 );
00054     matchVec.push_back( dummyMatch );
00055   }
00056 
00057   // declare auto_ptr for products
00058   std::auto_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >
00059     pOut( new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > );
00060   std::auto_ptr<int> pKey(new int);
00061 
00062   // go through given vector of jet combinations
00063   unsigned int idMatch = 0;
00064   typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
00065   for(MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
00066 
00067     // reset pointers
00068     resetCandidates();
00069 
00070     // build hypothesis
00071     buildHypo(evt, leps, mets, jets, *match, idMatch++);
00072 
00073     pOut->push_back( std::make_pair(hypo(), *match) );
00074   }
00075 
00076   // feed out hyps and matches
00077   evt.put(pOut);
00078 
00079   // build and feed out key
00080   buildKey();
00081   *pKey=key();
00082   evt.put(pKey, "Key");
00083 }

void TtSemiLepHypothesis::resetCandidates (  )  [protected]

reset candidate pointers before hypo build process

Definition at line 86 of file TtSemiLepHypothesis.cc.

References hadronicB_, lepton_, leptonicB_, lightQ_, lightQBar_, and neutrino_.

Referenced by produce().

00087 {
00088   lightQ_    = 0;
00089   lightQBar_ = 0;
00090   hadronicB_ = 0;
00091   leptonicB_ = 0;
00092   neutrino_  = 0;
00093   lepton_    = 0;
00094 }

template<typename C>
void TtSemiLepHypothesis::setCandidate ( const edm::Handle< C > &  handle,
const int idx,
reco::ShallowClonePtrCandidate *&  clone 
) [inline, protected]

use one object in a collection to set a ShallowClonePtrCandidate

Definition at line 80 of file TtSemiLepHypothesis.h.

References ptr.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), and TtSemiLepHypMVADisc::buildHypo().

00080                                                                                                                   {
00081   typedef typename C::value_type O;
00082   edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
00083   clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4(), ptr->vertex() );
00084 }


Member Data Documentation

bool TtSemiLepHypothesis::getMatch_ [protected]

Definition at line 59 of file TtSemiLepHypothesis.h.

Referenced by produce(), and TtSemiLepHypothesis().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::hadronicB_ [protected]

Definition at line 70 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().

edm::InputTag TtSemiLepHypothesis::jets_ [protected]

Definition at line 61 of file TtSemiLepHypothesis.h.

Referenced by produce().

int TtSemiLepHypothesis::key_ [protected]

Definition at line 66 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGeom::buildKey(), TtSemiLepHypGenMatch::buildKey(), TtSemiLepHypKinFit::buildKey(), TtSemiLepHypMVADisc::buildKey(), TtSemiLepHypWMassMaxSumPt::buildKey(), TtSemiLepHypMaxSumPtWMass::buildKey(), and key().

edm::InputTag TtSemiLepHypothesis::leps_ [protected]

Definition at line 62 of file TtSemiLepHypothesis.h.

Referenced by produce().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::lepton_ [protected]

Definition at line 73 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::leptonicB_ [protected]

Definition at line 71 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::lightQ_ [protected]

Definition at line 68 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::lightQBar_ [protected]

Definition at line 69 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().

edm::InputTag TtSemiLepHypothesis::match_ [protected]

Definition at line 64 of file TtSemiLepHypothesis.h.

Referenced by produce(), and TtSemiLepHypothesis().

edm::InputTag TtSemiLepHypothesis::mets_ [protected]

Definition at line 63 of file TtSemiLepHypothesis.h.

Referenced by produce().

reco::ShallowClonePtrCandidate* TtSemiLepHypothesis::neutrino_ [protected]

Definition at line 72 of file TtSemiLepHypothesis.h.

Referenced by TtSemiLepHypGenMatch::buildHypo(), TtSemiLepHypGeom::buildHypo(), TtSemiLepHypKinFit::buildHypo(), TtSemiLepHypWMassMaxSumPt::buildHypo(), TtSemiLepHypMaxSumPtWMass::buildHypo(), TtSemiLepHypMVADisc::buildHypo(), hypo(), resetCandidates(), and ~TtSemiLepHypothesis().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:45 2009 for CMSSW by  doxygen 1.5.4