#include <TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepHypWMassMaxSumPt.h>
Public Member Functions | |
TtSemiLepHypWMassMaxSumPt (const edm::ParameterSet &) | |
~TtSemiLepHypWMassMaxSumPt () | |
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 | |
int | maxNJets_ |
double | wMass_ |
Definition at line 7 of file TtSemiLepHypWMassMaxSumPt.h.
TtSemiLepHypWMassMaxSumPt::TtSemiLepHypWMassMaxSumPt | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 5 of file TtSemiLepHypWMassMaxSumPt.cc.
References Exception, and maxNJets_.
00005 : 00006 TtSemiLepHypothesis( cfg ), 00007 maxNJets_(cfg.getParameter<int> ("maxNJets")), 00008 wMass_ (cfg.getParameter<double>("wMass" )) 00009 { 00010 if(maxNJets_<4 && maxNJets_!=-1) 00011 throw cms::Exception("WrongConfig") 00012 << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n" 00013 << "It has to be larger than 4 or can be set to -1 to take all jets."; 00014 }
TtSemiLepHypWMassMaxSumPt::~TtSemiLepHypWMassMaxSumPt | ( | ) |
void TtSemiLepHypWMassMaxSumPt::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 19 of file TtSemiLepHypWMassMaxSumPt.cc.
References TtSemiLepEvtPartons::HadB, TtSemiLepHypothesis::hadronicB_, i, TtSemiLepHypothesis::isValid(), pfTauBenchmarkGeneric_cfi::jets, TtSemiLepEvtPartons::LepB, TtSemiLepEvtPartons::Lepton, TtSemiLepHypothesis::lepton_, TtSemiLepHypothesis::leptonicB_, TtSemiLepEvtPartons::LightQ, TtSemiLepHypothesis::lightQ_, TtSemiLepEvtPartons::LightQBar, TtSemiLepHypothesis::lightQBar_, maxNJets_, TtSemiLepHypothesis::neutrino_, TtSemiLepHypothesis::setCandidate(), sum(), and wMass_.
00024 { 00025 if(leps->empty() || mets->empty() || jets->size()<4){ 00026 // create empty hypothesis 00027 return; 00028 } 00029 00030 unsigned maxNJets = maxNJets_; 00031 if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size(); 00032 00033 match.clear(); 00034 for(unsigned int i=0; i<5; ++i) 00035 match.push_back(-1); 00036 00037 // ----------------------------------------------------- 00038 // associate those jets that get closest to the W mass 00039 // with their invariant mass to the hadronic W boson 00040 // ----------------------------------------------------- 00041 double wDist =-1.; 00042 std::vector<unsigned> closestToWMassIndices; 00043 for(unsigned idx=0; idx<maxNJets; ++idx){ 00044 for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){ 00045 reco::Particle::LorentzVector sum = 00046 (*jets)[idx].p4()+ 00047 (*jets)[jdx].p4(); 00048 if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){ 00049 wDist=fabs(sum.mass()-wMass_); 00050 closestToWMassIndices.clear(); 00051 closestToWMassIndices.push_back(idx); 00052 closestToWMassIndices.push_back(jdx); 00053 } 00054 } 00055 } 00056 00057 // ----------------------------------------------------- 00058 // associate those jets with maximum pt of the vectorial 00059 // sum to the hadronic decay chain 00060 // ----------------------------------------------------- 00061 double maxPt=-1.; 00062 unsigned hadB=0; 00063 for(unsigned idx=0; idx<maxNJets; ++idx){ 00064 // make sure it's not used up already from the hadronic W 00065 if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] ){ 00066 reco::Particle::LorentzVector sum = 00067 (*jets)[closestToWMassIndices[0]].p4()+ 00068 (*jets)[closestToWMassIndices[1]].p4()+ 00069 (*jets)[idx].p4(); 00070 if( maxPt<0. || maxPt<sum.pt() ){ 00071 maxPt=sum.pt(); 00072 hadB=idx; 00073 } 00074 } 00075 } 00076 00077 // ----------------------------------------------------- 00078 // associate the remaining jet with maximum pt of the 00079 // vectorial sum with the leading lepton with the 00080 // leptonic b quark 00081 // ----------------------------------------------------- 00082 maxPt=-1.; 00083 unsigned lepB=0; 00084 for(unsigned idx=0; idx<maxNJets; ++idx){ 00085 // make sure it's not used up already from the hadronic decay chain 00086 if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] && idx!=hadB) { 00087 reco::Particle::LorentzVector sum = 00088 (*jets)[idx].p4()+(*leps)[ 0 ].p4(); 00089 if( maxPt<0. || maxPt<sum.pt() ){ 00090 maxPt=sum.pt(); 00091 lepB=idx; 00092 } 00093 } 00094 } 00095 00096 // ----------------------------------------------------- 00097 // add jets 00098 // ----------------------------------------------------- 00099 if( isValid(closestToWMassIndices[0], jets) ){ 00100 setCandidate(jets, closestToWMassIndices[0], lightQ_); 00101 match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0]; 00102 } 00103 00104 if( isValid(closestToWMassIndices[1], jets) ){ 00105 setCandidate(jets, closestToWMassIndices[1], lightQBar_); 00106 match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1]; 00107 } 00108 00109 if( isValid(hadB, jets) ){ 00110 setCandidate(jets, hadB, hadronicB_); 00111 match[TtSemiLepEvtPartons::HadB] = hadB; 00112 } 00113 00114 if( isValid(lepB, jets) ){ 00115 setCandidate(jets, lepB, leptonicB_); 00116 match[TtSemiLepEvtPartons::LepB] = lepB; 00117 } 00118 00119 // ----------------------------------------------------- 00120 // add lepton 00121 // ----------------------------------------------------- 00122 setCandidate(leps, 0, lepton_); 00123 match[TtSemiLepEvtPartons::Lepton] = 0; 00124 00125 // ----------------------------------------------------- 00126 // add neutrino 00127 // ----------------------------------------------------- 00128 setCandidate(mets, 0, neutrino_); 00129 }
virtual void TtSemiLepHypWMassMaxSumPt::buildKey | ( | ) | [inline, private, virtual] |
build the event hypothesis key
Implements TtSemiLepHypothesis.
Definition at line 17 of file TtSemiLepHypWMassMaxSumPt.h.
References TtSemiLepHypothesis::key_, and TtEvent::kWMassMaxSumPt.
00017 { key_= TtSemiLeptonicEvent::kWMassMaxSumPt; };
int TtSemiLepHypWMassMaxSumPt::maxNJets_ [private] |
Definition at line 27 of file TtSemiLepHypWMassMaxSumPt.h.
Referenced by buildHypo(), and TtSemiLepHypWMassMaxSumPt().
double TtSemiLepHypWMassMaxSumPt::wMass_ [private] |