CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes
GenJetParticleSelector Class Reference

#include <GenJetParticleSelector.h>

Public Member Functions

 GenJetParticleSelector (const edm::ParameterSet &, edm::ConsumesCollector &iC)
 
void init (const edm::EventSetup &)
 
bool operator() (const reco::Candidate &)
 

Private Types

typedef std::vector< PdtEntryvpdt
 

Private Attributes

bool bInclude_
 
bool partons_
 
vpdt pdtList_
 
std::set< int > pIds_
 
bool stableOnly_
 
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordtableToken_
 

Detailed Description

Definition at line 25 of file GenJetParticleSelector.h.

Member Typedef Documentation

◆ vpdt

typedef std::vector<PdtEntry> GenJetParticleSelector::vpdt
private

Definition at line 32 of file GenJetParticleSelector.h.

Constructor & Destructor Documentation

◆ GenJetParticleSelector()

GenJetParticleSelector::GenJetParticleSelector ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 11 of file GenJetParticleSelector.cc.

References bInclude_, looper::cfg, Exception, genCandidates_cfi::excludeList, spr::find(), newFWLiteAna::found, partons_, pdtList_, and stableOnly_.

12  : stableOnly_(cfg.getParameter<bool>("stableOnly")),
13  partons_(false),
14  bInclude_(false),
15  tableToken_(iC.esConsumes()) {
16  const string excludeString("excludeList");
17  const string includeString("includeList");
18  vpdt includeList, excludeList;
19  vector<string> vPdtParams = cfg.getParameterNamesForType<vpdt>();
20  bool found = std::find(vPdtParams.begin(), vPdtParams.end(), includeString) != vPdtParams.end();
21  if (found)
22  includeList = cfg.getParameter<vpdt>(includeString);
23  found = find(vPdtParams.begin(), vPdtParams.end(), excludeString) != vPdtParams.end();
24  if (found)
25  excludeList = cfg.getParameter<vpdt>(excludeString);
26  const string partonsString("partons");
27  vector<string> vBoolParams = cfg.getParameterNamesForType<bool>();
28  found = find(vBoolParams.begin(), vBoolParams.end(), partonsString) != vBoolParams.end();
29  if (found)
30  partons_ = cfg.getParameter<bool>(partonsString);
31  bool bExclude = false;
32  if (!includeList.empty())
33  bInclude_ = true;
34  if (!excludeList.empty())
35  bExclude = true;
36 
37  if (bInclude_ && bExclude) {
38  throw cms::Exception("ConfigError", "not allowed to use both includeList and excludeList at the same time\n");
39  } else if (bInclude_) {
40  pdtList_ = includeList;
41  } else {
43  }
44  if (stableOnly_ && partons_) {
45  throw cms::Exception("ConfigError", "not allowed to have both stableOnly and partons true at the same time\n");
46  }
47 }
std::vector< PdtEntry > vpdt
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_

Member Function Documentation

◆ init()

void GenJetParticleSelector::init ( const edm::EventSetup es)

Definition at line 61 of file GenJetParticleSelector.cc.

References edm::EventSetup::getData(), mps_fire::i, pdtList_, and tableToken_.

Referenced by reco::modules::GenJetParticleSelectorEventSetupInit::init().

61  {
62  auto const& pdt = es.getData(tableToken_);
63 
64  for (vpdt::iterator i = pdtList_.begin(); i != pdtList_.end(); ++i)
65  i->setup(pdt);
66 }
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_

◆ operator()()

bool GenJetParticleSelector::operator() ( const reco::Candidate p)

Definition at line 49 of file GenJetParticleSelector.cc.

References funct::abs(), bInclude_, AlCaHLTBitMon_ParallelJobs::p, partons_, pIds_, stableOnly_, and mps_update::status.

49  {
50  int status = p.status();
51  int id = abs(p.pdgId());
52  if ((!stableOnly_ || status == 1) && !partons_ && ((pIds_.find(id) == pIds_.end()) ^ bInclude_))
53  return true;
54  else if (partons_ && (p.numberOfDaughters() > 0 && (p.daughter(0)->pdgId() == 91 || p.daughter(0)->pdgId() == 92)) &&
55  (((pIds_.find(id) == pIds_.end()) ^ bInclude_)))
56  return true;
57  else
58  return false;
59 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

Member Data Documentation

◆ bInclude_

bool GenJetParticleSelector::bInclude_
private

Definition at line 36 of file GenJetParticleSelector.h.

Referenced by GenJetParticleSelector(), and operator()().

◆ partons_

bool GenJetParticleSelector::partons_
private

Definition at line 34 of file GenJetParticleSelector.h.

Referenced by GenJetParticleSelector(), and operator()().

◆ pdtList_

vpdt GenJetParticleSelector::pdtList_
private

Definition at line 35 of file GenJetParticleSelector.h.

Referenced by GenJetParticleSelector(), and init().

◆ pIds_

std::set<int> GenJetParticleSelector::pIds_
private

Definition at line 37 of file GenJetParticleSelector.h.

Referenced by operator()().

◆ stableOnly_

bool GenJetParticleSelector::stableOnly_
private

Definition at line 33 of file GenJetParticleSelector.h.

Referenced by GenJetParticleSelector(), and operator()().

◆ tableToken_

const edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> GenJetParticleSelector::tableToken_
private

Definition at line 38 of file GenJetParticleSelector.h.

Referenced by init().