CMS 3D CMS Logo

GenJetParticleSelector.cc
Go to the documentation of this file.
6 #include <algorithm>
7 using namespace std;
8 using namespace edm;
9 
11  : stableOnly_(cfg.getParameter<bool>("stableOnly")), partons_(false), bInclude_(false) {
12  const string excludeString("excludeList");
13  const string includeString("includeList");
14  vpdt includeList, excludeList;
15  vector<string> vPdtParams = cfg.getParameterNamesForType<vpdt>();
16  bool found = std::find(vPdtParams.begin(), vPdtParams.end(), includeString) != vPdtParams.end();
17  if (found)
18  includeList = cfg.getParameter<vpdt>(includeString);
19  found = find(vPdtParams.begin(), vPdtParams.end(), excludeString) != vPdtParams.end();
20  if (found)
21  excludeList = cfg.getParameter<vpdt>(excludeString);
22  const string partonsString("partons");
23  vector<string> vBoolParams = cfg.getParameterNamesForType<bool>();
24  found = find(vBoolParams.begin(), vBoolParams.end(), partonsString) != vBoolParams.end();
25  if (found)
26  partons_ = cfg.getParameter<bool>(partonsString);
27  bool bExclude = false;
28  if (!includeList.empty())
29  bInclude_ = true;
30  if (!excludeList.empty())
31  bExclude = true;
32 
33  if (bInclude_ && bExclude) {
34  throw cms::Exception("ConfigError", "not allowed to use both includeList and excludeList at the same time\n");
35  } else if (bInclude_) {
36  pdtList_ = includeList;
37  } else {
39  }
40  if (stableOnly_ && partons_) {
41  throw cms::Exception("ConfigError", "not allowed to have both stableOnly and partons true at the same time\n");
42  }
43 }
44 
46  int status = p.status();
47  int id = abs(p.pdgId());
48  if ((!stableOnly_ || status == 1) && !partons_ && ((pIds_.find(id) == pIds_.end()) ^ bInclude_))
49  return true;
50  else if (partons_ && (p.numberOfDaughters() > 0 && (p.daughter(0)->pdgId() == 91 || p.daughter(0)->pdgId() == 92)) &&
51  (((pIds_.find(id) == pIds_.end()) ^ bInclude_)))
52  return true;
53  else
54  return false;
55 }
56 
58  for (vpdt::iterator i = pdtList_.begin(); i != pdtList_.end(); ++i)
59  i->setup(es);
60 }
GenJetParticleSelector::operator()
bool operator()(const reco::Candidate &)
Definition: GenJetParticleSelector.cc:45
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
funct::false
false
Definition: Factorize.h:34
genCandidates_cfi.excludeList
excludeList
Definition: genCandidates_cfi.py:7
mps_update.status
status
Definition: mps_update.py:69
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
newFWLiteAna.found
found
Definition: newFWLiteAna.py:118
EDMException.h
GenJetParticleSelector::partons_
bool partons_
Definition: GenJetParticleSelector.h:32
GenJetParticleSelector::pdtList_
vpdt pdtList_
Definition: GenJetParticleSelector.h:33
GenJetParticleSelector::bInclude_
bool bInclude_
Definition: GenJetParticleSelector.h:34
edm::ParameterSet
Definition: ParameterSet.h:36
PdtEntry.h
edm::EventSetup
Definition: EventSetup.h:57
GenJetParticleSelector::GenJetParticleSelector
GenJetParticleSelector(const edm::ParameterSet &, edm::ConsumesCollector &iC)
Definition: GenJetParticleSelector.cc:10
looper.cfg
cfg
Definition: looper.py:297
reco::Candidate
Definition: Candidate.h:27
GenJetParticleSelector::pIds_
std::set< int > pIds_
Definition: GenJetParticleSelector.h:35
std
Definition: JetResolutionObject.h:76
GenJetParticleSelector::init
void init(const edm::EventSetup &)
Definition: GenJetParticleSelector.cc:57
Exception
Definition: hltDiff.cc:246
GenJetParticleSelector.h
Candidate.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
GenJetParticleSelector::stableOnly_
bool stableOnly_
Definition: GenJetParticleSelector.h:31
GenJetParticleSelector::vpdt
std::vector< PdtEntry > vpdt
Definition: GenJetParticleSelector.h:30
edm::ConsumesCollector
Definition: ConsumesCollector.h:39