File indexing completed on 2023-03-17 10:45:17
0001 #include "CommonTools/CandAlgos/interface/GenJetParticleSelector.h"
0002 #include "DataFormats/Candidate/interface/Candidate.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "SimGeneral/HepPDTRecord/interface/PdtEntry.h"
0005 #include "FWCore/Utilities/interface/EDMException.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include <algorithm>
0008 using namespace std;
0009 using namespace edm;
0010
0011 GenJetParticleSelector::GenJetParticleSelector(const ParameterSet& cfg, edm::ConsumesCollector& iC)
0012 : stableOnly_(cfg.getParameter<bool>("stableOnly")),
0013 partons_(false),
0014 bInclude_(false),
0015 tableToken_(iC.esConsumes()) {
0016 const string excludeString("excludeList");
0017 const string includeString("includeList");
0018 vpdt includeList, excludeList;
0019 vector<string> vPdtParams = cfg.getParameterNamesForType<vpdt>();
0020 bool found = std::find(vPdtParams.begin(), vPdtParams.end(), includeString) != vPdtParams.end();
0021 if (found)
0022 includeList = cfg.getParameter<vpdt>(includeString);
0023 found = find(vPdtParams.begin(), vPdtParams.end(), excludeString) != vPdtParams.end();
0024 if (found)
0025 excludeList = cfg.getParameter<vpdt>(excludeString);
0026 const string partonsString("partons");
0027 vector<string> vBoolParams = cfg.getParameterNamesForType<bool>();
0028 found = find(vBoolParams.begin(), vBoolParams.end(), partonsString) != vBoolParams.end();
0029 if (found)
0030 partons_ = cfg.getParameter<bool>(partonsString);
0031 bool bExclude = false;
0032 if (!includeList.empty())
0033 bInclude_ = true;
0034 if (!excludeList.empty())
0035 bExclude = true;
0036
0037 if (bInclude_ && bExclude) {
0038 throw cms::Exception("ConfigError", "not allowed to use both includeList and excludeList at the same time\n");
0039 } else if (bInclude_) {
0040 pdtList_ = includeList;
0041 } else {
0042 pdtList_ = excludeList;
0043 }
0044 if (stableOnly_ && partons_) {
0045 throw cms::Exception("ConfigError", "not allowed to have both stableOnly and partons true at the same time\n");
0046 }
0047 }
0048
0049 bool GenJetParticleSelector::operator()(const reco::Candidate& p) {
0050 int status = p.status();
0051 int id = abs(p.pdgId());
0052 if ((!stableOnly_ || status == 1) && !partons_ && ((pIds_.find(id) == pIds_.end()) ^ bInclude_))
0053 return true;
0054 else if (partons_ && (p.numberOfDaughters() > 0 && (p.daughter(0)->pdgId() == 91 || p.daughter(0)->pdgId() == 92)) &&
0055 (((pIds_.find(id) == pIds_.end()) ^ bInclude_)))
0056 return true;
0057 else
0058 return false;
0059 }
0060
0061 void GenJetParticleSelector::init(const edm::EventSetup& es) {
0062 auto const& pdt = es.getData(tableToken_);
0063
0064 for (vpdt::iterator i = pdtList_.begin(); i != pdtList_.end(); ++i)
0065 i->setup(pdt);
0066 }