#include <GenJetParticleSelector.h>
Public Member Functions | |
GenJetParticleSelector (const edm::ParameterSet &) | |
void | init (const edm::EventSetup &) |
bool | operator() (const reco::Candidate &) |
Private Types | |
typedef std::vector< PdtEntry > | vpdt |
Private Attributes | |
bool | bInclude_ |
bool | partons_ |
vpdt | pdtList_ |
std::set< int > | pIds_ |
bool | stableOnly_ |
Definition at line 16 of file GenJetParticleSelector.h.
typedef std::vector<PdtEntry> GenJetParticleSelector::vpdt [private] |
Definition at line 22 of file GenJetParticleSelector.h.
GenJetParticleSelector::GenJetParticleSelector | ( | const edm::ParameterSet & | cfg | ) |
Definition at line 10 of file GenJetParticleSelector.cc.
References bInclude_, Exception, spr::find(), newFWLiteAna::found, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNamesForType(), partons_, pdtList_, and stableOnly_.
: stableOnly_(cfg.getParameter<bool>("stableOnly")), partons_(false), bInclude_(false) { const string excludeString("excludeList"); const string includeString("includeList"); vpdt includeList, excludeList; vector<string> vPdtParams = cfg.getParameterNamesForType<vpdt>(); bool found = std::find(vPdtParams.begin(), vPdtParams.end(), includeString) != vPdtParams.end(); if(found) includeList = cfg.getParameter<vpdt>(includeString); found = find(vPdtParams.begin(), vPdtParams.end(), excludeString) != vPdtParams.end(); if(found) excludeList = cfg.getParameter<vpdt>(excludeString); const string partonsString("partons"); vector<string> vBoolParams = cfg.getParameterNamesForType<bool>(); found = find(vBoolParams.begin(), vBoolParams.end(), partonsString) != vBoolParams.end(); if(found) partons_ = cfg.getParameter<bool>(partonsString); bool bExclude = false; if (includeList.size() > 0) bInclude_ = true; if (excludeList.size() > 0) bExclude = true; if (bInclude_ && bExclude) { throw cms::Exception("ConfigError", "not allowed to use both includeList and excludeList at the same time\n"); } else if (bInclude_) { pdtList_ = includeList; } else { pdtList_ = excludeList; } if(stableOnly_ && partons_) { throw cms::Exception("ConfigError", "not allowed to have both stableOnly and partons true at the same time\n"); } }
void GenJetParticleSelector::init | ( | const edm::EventSetup & | es | ) |
bool GenJetParticleSelector::operator() | ( | const reco::Candidate & | p | ) |
Definition at line 43 of file GenJetParticleSelector.cc.
References abs, bInclude_, reco::Candidate::daughter(), reco::Candidate::numberOfDaughters(), partons_, reco::Candidate::pdgId(), pIds_, stableOnly_, ntuplemaker::status, and reco::Candidate::status().
{ int status = p.status(); int id = abs(p.pdgId()); if((!stableOnly_ || status == 1) && !partons_ && ( (pIds_.find(id) == pIds_.end()) ^ bInclude_)) return true; else if(partons_ && (p.numberOfDaughters() > 0 && (p.daughter(0)->pdgId() == 91 || p.daughter(0)->pdgId() == 92)) && ( ((pIds_.find(id) == pIds_.end()) ^ bInclude_))) return true; else return false; }
bool GenJetParticleSelector::bInclude_ [private] |
Definition at line 26 of file GenJetParticleSelector.h.
Referenced by GenJetParticleSelector(), and operator()().
bool GenJetParticleSelector::partons_ [private] |
Definition at line 24 of file GenJetParticleSelector.h.
Referenced by GenJetParticleSelector(), and operator()().
vpdt GenJetParticleSelector::pdtList_ [private] |
Definition at line 25 of file GenJetParticleSelector.h.
Referenced by GenJetParticleSelector(), and init().
std::set<int> GenJetParticleSelector::pIds_ [private] |
Definition at line 27 of file GenJetParticleSelector.h.
Referenced by operator()().
bool GenJetParticleSelector::stableOnly_ [private] |
Definition at line 23 of file GenJetParticleSelector.h.
Referenced by GenJetParticleSelector(), and operator()().