CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
L1HPSPFTauQualityCut.cc
Go to the documentation of this file.
2 #include "FWCore/Utilities/interface/Exception.h" // cms::Exception
3 
5  : debug_(cfg.getUntrackedParameter<bool>("debug", false)) {
6  std::string pfCandTypeString = cfg.getParameter<std::string>("pfCandType");
7  if (pfCandTypeString == "chargedHadron")
9  else if (pfCandTypeString == "electron")
11  else if (pfCandTypeString == "muon")
13  else if (pfCandTypeString == "neutralHadron")
15  else if (pfCandTypeString == "photon")
17  else
18  throw cms::Exception("L1HPSPFTauQualityCut")
19  << "Invalid Configuration parameter 'pfCandType' = '" << pfCandTypeString << "' !!\n";
20 
21  std::string dzCutString = cfg.getParameter<std::string>("dzCut");
22  if (dzCutString == "disabled")
23  dzCut_ = kDisabled;
24  else if (dzCutString == "enabled_primary")
26  else if (dzCutString == "enabled_pileup")
28  else
29  throw cms::Exception("L1HPSPFTauQualityCut")
30  << "Invalid Configuration parameter 'dzCut' = '" << dzCutString << "' !!\n";
31 
32  minPt_ = cfg.getParameter<double>("minPt");
33  maxDz_ = (cfg.exists("maxDz")) ? cfg.getParameter<double>("maxDz") : 1.e+3;
34 }
35 
36 bool L1HPSPFTauQualityCut::operator()(const l1t::PFCandidate& pfCand, float_t primaryVertex_z) const {
37  if (pfCand.id() == pfCandType_) {
38  if (pfCand.pt() < minPt_) {
39  return false;
40  }
41 
42  if (pfCand.charge() != 0) {
44  const l1t::PFTrackRef& pfCand_track = pfCand.pfTrack();
45  double dz = std::fabs(pfCand_track->vertex().z() - primaryVertex_z);
46  if (dzCut_ == kEnabledPrimary && dz > maxDz_)
47  return false;
48  if (dzCut_ == kEnabledPileup && dz <= maxDz_)
49  return false;
50  }
51  } else if (dzCut_ == kEnabledPileup) {
52  return false; // CV: only consider charged PFCands as originating from pileup
53  }
54  }
55  return true;
56 }
57 
59 
60 int L1HPSPFTauQualityCut::dzCut() const { return dzCut_; }
61 
62 float_t L1HPSPFTauQualityCut::minPt() const { return minPt_; }
63 
64 float_t L1HPSPFTauQualityCut::maxDz() const { return maxDz_; }
65 
67  const std::string& pfCandType,
68  const std::string& dzCut,
69  bool debug) {
70  edm::ParameterSet cfg_pfCandType = cfg.getParameter<edm::ParameterSet>(pfCandType);
71  cfg_pfCandType.addParameter<std::string>("pfCandType", pfCandType);
72  cfg_pfCandType.addParameter<std::string>("dzCut", dzCut);
73  cfg_pfCandType.addUntrackedParameter<bool>("debug", debug);
74  L1HPSPFTauQualityCut qualityCut(cfg_pfCandType);
75  return qualityCut;
76 }
77 
78 std::vector<L1HPSPFTauQualityCut> readL1PFTauQualityCuts(const edm::ParameterSet& cfg,
79  const std::string& dzCut,
80  bool debug) {
81  std::vector<L1HPSPFTauQualityCut> qualityCuts;
82  qualityCuts.push_back(readL1PFTauQualityCut(cfg, "chargedHadron", dzCut, debug));
83  qualityCuts.push_back(readL1PFTauQualityCut(cfg, "electron", dzCut, debug));
84  qualityCuts.push_back(readL1PFTauQualityCut(cfg, "muon", dzCut, debug));
85  qualityCuts.push_back(readL1PFTauQualityCut(cfg, "photon", dzCut, debug));
86  qualityCuts.push_back(readL1PFTauQualityCut(cfg, "neutralHadron", dzCut, debug));
87  return qualityCuts;
88 }
89 
90 bool isSelected(const std::vector<L1HPSPFTauQualityCut>& qualityCuts,
91  const l1t::PFCandidate& pfCand,
92  float_t primaryVertex_z) {
93  for (auto qualityCut : qualityCuts) {
94  if (!qualityCut(pfCand, primaryVertex_z))
95  return false;
96  }
97  return true;
98 }
l1t::PFCandidate::ParticleType pfCandType_
tuple cfg
Definition: looper.py:296
double pt() const final
transverse momentum
bool operator()(const l1t::PFCandidate &pfCand, float_t primaryVertexZ) const
returns true (false) if PFCandidate passes (fails) quality cuts
bool exists(std::string const &parameterName) const
checks if a parameter exists
ParticleType id() const
Definition: PFCandidate.h:34
l1t::PFCandidate::ParticleType pfCandType() const
accessor functions
const PFTrackRef & pfTrack() const
Definition: PFCandidate.h:36
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
std::vector< L1HPSPFTauQualityCut > readL1PFTauQualityCuts(const edm::ParameterSet &cfg, const std::string &dzCut, bool debug=false)
#define debug
Definition: HDRShower.cc:19
void addUntrackedParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:192
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool isSelected(const std::vector< L1HPSPFTauQualityCut > &qualityCuts, const l1t::PFCandidate &pfCand, float_t primaryVertexZ)
L1HPSPFTauQualityCut(const edm::ParameterSet &cfg)
constructor
L1HPSPFTauQualityCut readL1PFTauQualityCut(const edm::ParameterSet &cfg, const std::string &pfCandType, const std::string &dzCut, bool debug)
int charge() const final
electric charge