CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GenJetParticleSelector.cc
Go to the documentation of this file.
7 #include <algorithm>
8 using namespace std;
9 using namespace edm;
10 
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 {
42  pdtList_ = excludeList;
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 }
48 
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 }
60 
62  auto const& pdt = es.getData(tableToken_);
63 
64  for (vpdt::iterator i = pdtList_.begin(); i != pdtList_.end(); ++i)
65  i->setup(pdt);
66 }
GenJetParticleSelector(const edm::ParameterSet &, edm::ConsumesCollector &iC)
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
tuple cfg
Definition: looper.py:296
virtual int status() const =0
status word
std::vector< PdtEntry > vpdt
list status
Definition: mps_update.py:107
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
Definition: ParameterSet.h:179
bool getData(T &iHolder) const
Definition: EventSetup.h:122
virtual size_type numberOfDaughters() const =0
number of daughters
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void init(const edm::EventSetup &)
virtual int pdgId() const =0
PDG identifier.
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool operator()(const reco::Candidate &)
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283