#include <TtSemiLepJetCombWMassMaxSumPt.h>
Public Member Functions | |
TtSemiLepJetCombWMassMaxSumPt (const edm::ParameterSet &) | |
~TtSemiLepJetCombWMassMaxSumPt () | |
Private Member Functions | |
virtual void | beginJob () |
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_ |
double | wMass_ |
Definition at line 9 of file TtSemiLepJetCombWMassMaxSumPt.h.
TtSemiLepJetCombWMassMaxSumPt::TtSemiLepJetCombWMassMaxSumPt | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 6 of file TtSemiLepJetCombWMassMaxSumPt.cc.
References Exception, and maxNJets_.
: jets_ (cfg.getParameter<edm::InputTag>("jets" )), leps_ (cfg.getParameter<edm::InputTag>("leps" )), maxNJets_ (cfg.getParameter<int> ("maxNJets" )), wMass_ (cfg.getParameter<double> ("wMass" )), 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> > >(); produces<int>("NumberOfConsideredJets"); }
TtSemiLepJetCombWMassMaxSumPt::~TtSemiLepJetCombWMassMaxSumPt | ( | ) |
Definition at line 25 of file TtSemiLepJetCombWMassMaxSumPt.cc.
{ }
virtual void TtSemiLepJetCombWMassMaxSumPt::beginJob | ( | void | ) | [inline, private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 18 of file TtSemiLepJetCombWMassMaxSumPt.h.
{};
virtual void TtSemiLepJetCombWMassMaxSumPt::endJob | ( | void | ) | [inline, private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 20 of file TtSemiLepJetCombWMassMaxSumPt.h.
{};
bool TtSemiLepJetCombWMassMaxSumPt::isValid | ( | const int & | idx, |
const edm::Handle< std::vector< pat::Jet > > & | jets | ||
) | [inline, private] |
Definition at line 22 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
{ return (0<=idx && idx<(int)jets->size()); };
void TtSemiLepJetCombWMassMaxSumPt::produce | ( | edm::Event & | evt, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 30 of file TtSemiLepJetCombWMassMaxSumPt.cc.
References bTagAlgorithm_, edm::Event::getByLabel(), TtSemiLepDaughter::HadB, i, isValid(), fwrapper::jets, jets_, TtSemiLepDaughter::LepB, leps_, TtFullHadDaughter::LightQ, TtFullHadDaughter::LightQBar, match(), maxBDiscLightJets_, maxNJets_, minBDiscBJets_, edm::Event::put(), useBTagging_, and wMass_.
{ std::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >); std::auto_ptr<int> pJetsConsidered(new 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 without lepton candidate or less than 4 jets or no MET if(leps->empty() || jets->size() < 4){ pOut->push_back( match ); evt.put(pOut); *pJetsConsidered = jets->size(); evt.put(pJetsConsidered, "NumberOfConsideredJets"); return; } unsigned maxNJets = maxNJets_; if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size(); *pJetsConsidered = maxNJets; evt.put(pJetsConsidered, "NumberOfConsideredJets"); 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 jets that get closest to the W mass // with their invariant mass to the hadronic W boson // ----------------------------------------------------- double wDist =-1.; std::vector<int> closestToWMassIndices; closestToWMassIndices.push_back(-1); closestToWMassIndices.push_back(-1); for(unsigned idx=0; idx<maxNJets; ++idx){ if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue; for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){ if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue; reco::Particle::LorentzVector sum = (*jets)[idx].p4()+ (*jets)[jdx].p4(); if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){ wDist=fabs(sum.mass()-wMass_); closestToWMassIndices.clear(); closestToWMassIndices.push_back(idx); closestToWMassIndices.push_back(jdx); } } } // ----------------------------------------------------- // associate those jets with maximum pt of the vectorial // sum to the hadronic decay chain // ----------------------------------------------------- double maxPt=-1.; int hadB=-1; if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) { for(unsigned 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!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] ){ reco::Particle::LorentzVector sum = (*jets)[closestToWMassIndices[0]].p4()+ (*jets)[closestToWMassIndices[1]].p4()+ (*jets)[idx].p4(); if( maxPt<0. || maxPt<sum.pt() ){ maxPt=sum.pt(); hadB=idx; } } } } // ----------------------------------------------------- // associate the remaining jet with maximum pt of the // vectorial sum with the leading lepton with the // leptonic b quark // ----------------------------------------------------- maxPt=-1.; int lepB=-1; for(unsigned 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!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] && (int)idx!=hadB) { reco::Particle::LorentzVector sum = (*jets)[idx].p4()+(*leps)[ 0 ].p4(); if( maxPt<0. || maxPt<sum.pt() ){ maxPt=sum.pt(); lepB=idx; } } } match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0]; match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1]; match[TtSemiLepEvtPartons::HadB ] = hadB; match[TtSemiLepEvtPartons::LepB ] = lepB; pOut->push_back( match ); evt.put(pOut); }
std::string TtSemiLepJetCombWMassMaxSumPt::bTagAlgorithm_ [private] |
Definition at line 29 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
Definition at line 22 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
Definition at line 25 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
double TtSemiLepJetCombWMassMaxSumPt::maxBDiscLightJets_ [private] |
Definition at line 31 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
int TtSemiLepJetCombWMassMaxSumPt::maxNJets_ [private] |
Definition at line 26 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce(), and TtSemiLepJetCombWMassMaxSumPt().
double TtSemiLepJetCombWMassMaxSumPt::minBDiscBJets_ [private] |
Definition at line 30 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
bool TtSemiLepJetCombWMassMaxSumPt::useBTagging_ [private] |
Definition at line 28 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().
double TtSemiLepJetCombWMassMaxSumPt::wMass_ [private] |
Definition at line 27 of file TtSemiLepJetCombWMassMaxSumPt.h.
Referenced by produce().