#include <TtSemiLepJetCombGeom.h>
Public Member Functions | |
TtSemiLepJetCombGeom (const edm::ParameterSet &) | |
~TtSemiLepJetCombGeom () | |
Private Member Functions | |
virtual void | beginJob () |
double | distance (const math::XYZTLorentzVector &, const math::XYZTLorentzVector &) |
virtual void | endJob () |
bool | isValid (const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets) |
virtual void | produce (edm::Event &evt, const edm::EventSetup &setup) |
Private Attributes | |
std::string | bTagAlgorithm_ |
edm::InputTag | jets_ |
edm::InputTag | leps_ |
double | maxBDiscLightJets_ |
int | maxNJets_ |
double | minBDiscBJets_ |
bool | useBTagging_ |
bool | useDeltaR_ |
Definition at line 10 of file TtSemiLepJetCombGeom.h.
TtSemiLepJetCombGeom::TtSemiLepJetCombGeom | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 9 of file TtSemiLepJetCombGeom.cc.
References Exception, and maxNJets_.
: jets_ (cfg.getParameter<edm::InputTag>("jets" )), leps_ (cfg.getParameter<edm::InputTag>("leps" )), maxNJets_ (cfg.getParameter<int> ("maxNJets" )), useDeltaR_ (cfg.getParameter<bool> ("useDeltaR" )), useBTagging_ (cfg.getParameter<bool> ("useBTagging" )), bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )), minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )), maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets")) { if(maxNJets_<4 && maxNJets_!=-1) throw cms::Exception("WrongConfig") << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n" << "It has to be larger than 4 or can be set to -1 to take all jets."; produces< std::vector<std::vector<int> > >(); }
TtSemiLepJetCombGeom::~TtSemiLepJetCombGeom | ( | ) |
Definition at line 27 of file TtSemiLepJetCombGeom.cc.
{ }
virtual void TtSemiLepJetCombGeom::beginJob | ( | void | ) | [inline, private, virtual] |
double TtSemiLepJetCombGeom::distance | ( | const math::XYZTLorentzVector & | v1, |
const math::XYZTLorentzVector & | v2 | ||
) | [private] |
Definition at line 141 of file TtSemiLepJetCombGeom.cc.
References useDeltaR_.
Referenced by produce().
{ // calculate the distance between two lorentz vectors // using DeltaR or DeltaTheta if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2); return fabs(v1.theta() - v2.theta()); }
virtual void TtSemiLepJetCombGeom::endJob | ( | void | ) | [inline, private, virtual] |
bool TtSemiLepJetCombGeom::isValid | ( | const int & | idx, |
const edm::Handle< std::vector< pat::Jet > > & | jets | ||
) | [inline, private] |
Definition at line 23 of file TtSemiLepJetCombGeom.h.
References analyzePatCleaning_cfg::jets.
Referenced by produce().
{ return (0<=idx && idx<(int)jets->size()); };
void TtSemiLepJetCombGeom::produce | ( | edm::Event & | evt, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 32 of file TtSemiLepJetCombGeom.cc.
References bTagAlgorithm_, distance(), edm::Event::getByLabel(), TtSemiLepDaughter::HadB, i, isValid(), analyzePatCleaning_cfg::jets, jets_, TtSemiLepDaughter::LepB, leps_, TtFullHadDaughter::LightQ, TtFullHadDaughter::LightQBar, match(), maxBDiscLightJets_, maxNJets_, minBDiscBJets_, p4, edm::Event::put(), and useBTagging_.
{ std::auto_ptr< std::vector<std::vector<int> > >pOut (new std::vector<std::vector<int> >); std::vector<int> match; for(unsigned int i = 0; i < 4; ++i) match.push_back( -1 ); // get jets edm::Handle< std::vector<pat::Jet> > jets; evt.getByLabel(jets_, jets); // get leptons edm::Handle< edm::View<reco::RecoCandidate> > leps; evt.getByLabel(leps_, leps); // skip events with without lepton candidate or less than 4 jets if(leps->empty() || jets->size() < 4){ pOut->push_back( match ); evt.put(pOut); return; } unsigned maxNJets = maxNJets_; if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size(); std::vector<bool> isBJet; std::vector<bool> isLJet; int cntBJets = 0; if(useBTagging_) { for(unsigned int idx=0; idx<maxNJets; ++idx) { isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) ); isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) ); if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++; } } // ----------------------------------------------------- // associate those two jets to the hadronic W boson that // have the smallest distance to each other // ----------------------------------------------------- double minDist=-1.; int lightQ =-1; int lightQBar=-1; for(unsigned int idx=0; idx<maxNJets; ++idx){ if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue; for(unsigned int jdx=(idx+1); jdx<maxNJets; ++jdx){ if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue; double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4()); if( minDist<0. || dist<minDist ){ minDist=dist; lightQ =idx; lightQBar=jdx; } } } reco::Particle::LorentzVector wHad; if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4(); // ----------------------------------------------------- // associate to the hadronic b quark the remaining jet // that has the smallest distance to the hadronic W // ----------------------------------------------------- minDist=-1.; int hadB=-1; if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) { for(unsigned int idx=0; idx<maxNJets; ++idx){ if(useBTagging_ && !isBJet[idx]) continue; // make sure it's not used up already from the hadronic W if( (int)idx!=lightQ && (int)idx!=lightQBar ) { double dist = distance((*jets)[idx].p4(), wHad); if( minDist<0. || dist<minDist ){ minDist=dist; hadB=idx; } } } } // ----------------------------------------------------- // associate to the leptonic b quark the remaining jet // that has the smallest distance to the leading lepton // ----------------------------------------------------- minDist=-1.; int lepB=-1; for(unsigned int idx=0; idx<maxNJets; ++idx){ if(useBTagging_ && !isBJet[idx]) continue; // make sure it's not used up already from the hadronic decay chain if( (int)idx!=lightQ && (int)idx!=lightQBar && (int)idx!=hadB ){ double dist = distance((*jets)[idx].p4(), (*leps)[0].p4()); if( minDist<0. || dist<minDist ){ minDist=dist; lepB=idx; } } } match[TtSemiLepEvtPartons::LightQ ] = lightQ; match[TtSemiLepEvtPartons::LightQBar] = lightQBar; match[TtSemiLepEvtPartons::HadB ] = hadB; match[TtSemiLepEvtPartons::LepB ] = lepB; pOut->push_back( match ); evt.put(pOut); }
std::string TtSemiLepJetCombGeom::bTagAlgorithm_ [private] |
Definition at line 31 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
edm::InputTag TtSemiLepJetCombGeom::jets_ [private] |
Definition at line 26 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
edm::InputTag TtSemiLepJetCombGeom::leps_ [private] |
Definition at line 27 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
double TtSemiLepJetCombGeom::maxBDiscLightJets_ [private] |
Definition at line 33 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
int TtSemiLepJetCombGeom::maxNJets_ [private] |
Definition at line 28 of file TtSemiLepJetCombGeom.h.
Referenced by produce(), and TtSemiLepJetCombGeom().
double TtSemiLepJetCombGeom::minBDiscBJets_ [private] |
Definition at line 32 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
bool TtSemiLepJetCombGeom::useBTagging_ [private] |
Definition at line 30 of file TtSemiLepJetCombGeom.h.
Referenced by produce().
bool TtSemiLepJetCombGeom::useDeltaR_ [private] |
Definition at line 29 of file TtSemiLepJetCombGeom.h.
Referenced by distance().