CMS 3D CMS Logo

L1GTSingleCollectionCut.h
Go to the documentation of this file.
1 #ifndef L1Trigger_Phase2L1GT_L1GTSingleCollectionCut_h
2 #define L1Trigger_Phase2L1GT_L1GTSingleCollectionCut_h
3 
8 
10 
11 #include "L1GTOptionalParam.h"
12 
13 #include <functional>
14 #include <optional>
15 #include <cstdint>
16 
17 namespace l1t {
18 
19  template <typename T, typename K>
20  inline std::vector<T> getParamVector(const std::string& name,
22  std::function<T(K)> conv) {
23  const std::vector<K>& values = config.getParameter<std::vector<K>>(name);
24  std::vector<T> convertedValues(values.size());
25  for (std::size_t i = 0; i < values.size(); i++) {
26  convertedValues[i] = conv(values[i]);
27  }
28  return convertedValues;
29  }
30 
32  public:
34  const edm::ParameterSet& lutConfig,
35  const L1GTScales& scales)
36  : scales_(scales),
37  tag_(config.getParameter<edm::InputTag>("tag")),
38  minPt_(getOptionalParam<int, double>(
39  "minPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
40  maxPt_(getOptionalParam<int, double>(
41  "maxPt", config, [&scales](double value) { return scales.to_hw_pT_ceil(value); })),
42  minEta_(getOptionalParam<int, double>(
43  "minEta", config, [&scales](double value) { return scales.to_hw_eta_floor(value); })),
44  maxEta_(getOptionalParam<int, double>(
45  "maxEta", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
46  minPhi_(getOptionalParam<int, double>(
47  "minPhi", config, [&scales](double value) { return scales.to_hw_phi_floor(value); })),
48  maxPhi_(getOptionalParam<int, double>(
49  "maxPhi", config, [&scales](double value) { return scales.to_hw_phi_ceil(value); })),
50  minZ0_(getOptionalParam<int, double>(
51  "minZ0", config, [&scales](double value) { return scales.to_hw_z0_floor(value); })),
52  maxZ0_(getOptionalParam<int, double>(
53  "maxZ0", config, [&scales](double value) { return scales.to_hw_z0_ceil(value); })),
54  minScalarSumPt_(getOptionalParam<int, double>(
55  "minScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
56  maxScalarSumPt_(getOptionalParam<int, double>(
57  "maxScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT_ceil(value); })),
58  minQualityScore_(getOptionalParam<unsigned int>("minQualityScore", config)),
59  maxQualityScore_(getOptionalParam<unsigned int>("maxQualityScore", config)),
60  qualityFlags_(getOptionalParam<unsigned int>("qualityFlags", config)),
61  minAbsEta_(getOptionalParam<int, double>(
62  "minAbsEta", config, [&scales](double value) { return scales.to_hw_eta_floor(value); })),
63  maxAbsEta_(getOptionalParam<int, double>(
64  "maxAbsEta", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
65  minIsolationPt_(getOptionalParam<int, double>(
66  "minRelIsolationPt", config, [&scales](double value) { return scales.to_hw_isolationPT_floor(value); })),
67  maxIsolationPt_(getOptionalParam<int, double>(
68  "maxRelIsolationPt", config, [&scales](double value) { return scales.to_hw_isolationPT_ceil(value); })),
69  minRelIsolationPt_(getOptionalParam<int, double>(
70  "minRelIsolationPt",
71  config,
72  [&scales](double value) { return scales.to_hw_relative_isolationPT_floor(value); })),
73  maxRelIsolationPt_(getOptionalParam<int, double>(
74  "maxRelIsolationPt",
75  config,
76  [&scales](double value) { return scales.to_hw_relative_isolationPT_ceil(value); })),
77  regionsAbsEtaLowerBounds_(getParamVector<int, double>(
78  "regionsAbsEtaLowerBounds", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
79  regionsMinPt_(getParamVector<int, double>(
80  "regionsMinPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
81  regionsMaxRelIsolationPt_(getParamVector<int, double>(
82  "regionsMaxRelIsolationPt",
83  config,
84  [&scales](double value) { return scales.to_hw_relative_isolationPT_ceil(value); })),
85  regionsQualityFlags_(config.getParameter<std::vector<unsigned int>>("regionsQualityFlags")),
86  minPrimVertDz_(getOptionalParam<int, double>(
87  "minPrimVertDz", config, [&scales](double value) { return scales.to_hw_z0_floor(value); })),
88  maxPrimVertDz_(getOptionalParam<int, double>(
89  "maxPrimVertDz", config, [&scales](double value) { return scales.to_hw_z0_ceil(value); })),
90  primVertex_(getOptionalParam<unsigned int>("primVertex", config)) {}
91 
92  bool checkObject(const P2GTCandidate& obj) const {
93  bool result = true;
94 
95  result &= minPt_ ? (obj.hwPT() > minPt_) : true;
96  result &= maxPt_ ? (obj.hwPT() < maxPt_) : true;
97 
98  result &= minEta_ ? (obj.hwEta() > minEta_) : true;
99  result &= maxEta_ ? (obj.hwEta() < maxEta_) : true;
100 
101  result &= minPhi_ ? (obj.hwPhi() > minPhi_) : true;
102  result &= maxPhi_ ? (obj.hwPhi() < maxPhi_) : true;
103 
104  result &= minZ0_ ? (obj.hwZ0() > minZ0_) : true;
105  result &= maxZ0_ ? (obj.hwZ0() < maxZ0_) : true;
106 
107  result &= minAbsEta_ ? (abs(obj.hwEta()) > minAbsEta_) : true;
108  result &= maxAbsEta_ ? (abs(obj.hwEta()) < maxAbsEta_) : true;
109 
110  result &= minScalarSumPt_ ? (obj.hwScalarSumPT() > minScalarSumPt_) : true;
111  result &= maxScalarSumPt_ ? (obj.hwScalarSumPT() < maxScalarSumPt_) : true;
112 
113  result &= minQualityScore_ ? (obj.hwQualityScore() > minQualityScore_) : true;
114  result &= maxQualityScore_ ? (obj.hwQualityScore() < maxQualityScore_) : true;
115  result &= qualityFlags_ ? (obj.hwQualityFlags().to_uint() & qualityFlags_.value()) == qualityFlags_ : true;
116 
117  result &= minIsolationPt_ ? (obj.hwIsolationPT() > minIsolationPt_) : true;
118  result &= maxIsolationPt_ ? (obj.hwIsolationPT() < maxIsolationPt_) : true;
119 
120  result &= minRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION >
121  static_cast<int64_t>(minRelIsolationPt_.value()) * obj.hwPT().to_int64()
122  : true;
123  result &= maxRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
124  static_cast<int64_t>(maxRelIsolationPt_.value()) * obj.hwPT().to_int64()
125  : true;
126 
128  return result;
129  }
130 
131  bool checkPrimaryVertices(const P2GTCandidate& obj, const P2GTCandidateCollection& primVertCol) const {
132  if (!minPrimVertDz_ && !maxPrimVertDz_) {
133  return true;
134  } else {
135  if (!primVertex_) {
136  throw cms::Exception("Configuration")
137  << "When a min/max primary vertex cut is set the primary vertex to cut\n"
138  << "on has to be specified with with config parameter \'primVertex\'. ";
139  }
140  if (primVertex_ < primVertCol.size()) {
141  uint32_t dZ = abs(obj.hwZ0() - primVertCol[primVertex_.value()].hwZ0());
142 
143  bool result = true;
144  result &= minPrimVertDz_ ? dZ > minPrimVertDz_ : true;
145  result &= maxPrimVertDz_ ? dZ < maxPrimVertDz_ : true;
146  return result;
147  } else {
148  return false;
149  }
150  }
151  }
152 
154  desc.add<edm::InputTag>("tag");
155  desc.addOptional<double>("minPt");
156  desc.addOptional<double>("maxPt");
157  desc.addOptional<double>("minEta");
158  desc.addOptional<double>("maxEta");
159  desc.addOptional<double>("minPhi");
160  desc.addOptional<double>("maxPhi");
161  desc.addOptional<double>("minZ0");
162  desc.addOptional<double>("maxZ0");
163  desc.addOptional<double>("minScalarSumPt");
164  desc.addOptional<double>("maxScalarSumPt");
165  desc.addOptional<unsigned int>("minQualityScore");
166  desc.addOptional<unsigned int>("maxQualityScore");
167  desc.addOptional<unsigned int>("qualityFlags");
168  desc.add<std::vector<unsigned int>>("regions", {});
169  desc.addOptional<double>("minAbsEta");
170  desc.addOptional<double>("maxAbsEta");
171  desc.addOptional<double>("minIsolationPt");
172  desc.addOptional<double>("maxIsolationPt");
173  desc.addOptional<double>("minRelIsolationPt");
174  desc.addOptional<double>("maxRelIsolationPt");
175  desc.add<std::vector<double>>("regionsAbsEtaLowerBounds", {});
176  desc.add<std::vector<double>>("regionsMinPt", {});
177  desc.add<std::vector<double>>("regionsMaxRelIsolationPt", {});
178  desc.add<std::vector<unsigned int>>("regionsQualityFlags", {});
179  desc.addOptional<double>("minPrimVertDz");
180  desc.addOptional<double>("maxPrimVertDz");
181  desc.addOptional<unsigned int>("primVertex");
182  }
183 
184  const edm::InputTag& tag() const { return tag_; }
185 
186  private:
188  bool res = true;
189 
190  unsigned int index;
191  index = atIndex(obj.hwEta());
192  res &= regionsMinPt_.empty() ? true : obj.hwPT() > regionsMinPt_[index];
193  res &= regionsMaxRelIsolationPt_.empty()
194  ? true
195  : obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
196  static_cast<int64_t>(regionsMaxRelIsolationPt_[index]) * obj.hwPT().to_int64();
197  res &= regionsQualityFlags_.empty()
198  ? true
199  : (obj.hwQualityFlags().to_uint() & regionsQualityFlags_[index]) == regionsQualityFlags_[index];
200  return res;
201  }
202 
203  unsigned int atIndex(int objeta) const {
204  // Function that checks at which index the eta of the object is
205  // If object abs(eta) < regionsAbsEtaLowerBounds_[0] the function returns the last index,
206  // Might be undesired behaviour
207  for (unsigned int i = 0; i < regionsAbsEtaLowerBounds_.size() - 1; i++) {
208  if (std::abs(objeta) >= regionsAbsEtaLowerBounds_[i] && std::abs(objeta) < regionsAbsEtaLowerBounds_[i + 1]) {
209  return i;
210  }
211  }
212  return regionsAbsEtaLowerBounds_.size() - 1;
213  }
214 
215  private:
218  const std::optional<int> minPt_;
219  const std::optional<int> maxPt_;
220  const std::optional<int> minEta_;
221  const std::optional<int> maxEta_;
222  const std::optional<int> minPhi_;
223  const std::optional<int> maxPhi_;
224  const std::optional<int> minZ0_;
225  const std::optional<int> maxZ0_;
226  const std::optional<int> minScalarSumPt_;
227  const std::optional<int> maxScalarSumPt_;
228  const std::optional<unsigned int> minQualityScore_;
229  const std::optional<unsigned int> maxQualityScore_;
230  const std::optional<unsigned int> qualityFlags_;
231  const std::optional<int> minAbsEta_;
232  const std::optional<int> maxAbsEta_;
233  const std::optional<int> minIsolationPt_;
234  const std::optional<int> maxIsolationPt_;
235  const std::optional<int> minRelIsolationPt_;
236  const std::optional<int> maxRelIsolationPt_;
237  const std::vector<int> regionsAbsEtaLowerBounds_;
238  const std::vector<int> regionsMinPt_;
239  const std::vector<int> regionsMaxRelIsolationPt_;
240  const std::vector<unsigned int> regionsQualityFlags_;
241  const std::optional<int> minPrimVertDz_;
242  const std::optional<int> maxPrimVertDz_;
243  const std::optional<unsigned int> primVertex_;
244  };
245 
246 } // namespace l1t
247 
248 #endif // L1Trigger_Phase2L1GT_L1GTSingleCollectionCut_h
const std::optional< int > maxAbsEta_
const std::optional< unsigned int > maxQualityScore_
const std::optional< int > maxRelIsolationPt_
bool checkEtaDependentCuts(const P2GTCandidate &obj) const
const std::vector< int > regionsMaxRelIsolationPt_
const std::optional< int > minPhi_
bool checkObject(const P2GTCandidate &obj) const
const std::optional< unsigned int > primVertex_
const std::optional< int > minPrimVertDz_
const std::vector< unsigned int > regionsQualityFlags_
std::vector< P2GTCandidate > P2GTCandidateCollection
Definition: P2GTCandidate.h:16
const std::optional< int > maxEta_
const std::optional< int > maxPhi_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
const edm::InputTag & tag() const
const std::vector< int > regionsAbsEtaLowerBounds_
const std::optional< int > maxScalarSumPt_
bool checkPrimaryVertices(const P2GTCandidate &obj, const P2GTCandidateCollection &primVertCol) const
delete x;
Definition: CaloConfig.h:22
Definition: config.py:1
const std::optional< unsigned int > qualityFlags_
L1GTSingleCollectionCut(const edm::ParameterSet &config, const edm::ParameterSet &lutConfig, const L1GTScales &scales)
Definition: Electron.h:6
const std::optional< int > minAbsEta_
std::optional< T > getOptionalParam(const std::string &name, const edm::ParameterSet &config, std::function< T(K)> conv)
unsigned int atIndex(int objeta) const
const std::optional< int > minRelIsolationPt_
const std::optional< int > minEta_
const std::optional< int > maxPt_
const std::optional< int > maxPrimVertDz_
const std::optional< int > minZ0_
const std::optional< int > maxIsolationPt_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: value.py:1
const std::optional< unsigned int > minQualityScore_
EPOS::IO_EPOS conv
static constexpr int RELATIVE_ISOLATION_RESOLUTION
Definition: L1GTScales.h:12
const std::vector< int > regionsMinPt_
HLT enums.
std::vector< T > getParamVector(const std::string &name, const edm::ParameterSet &config, std::function< T(K)> conv)
long double T
const std::optional< int > minIsolationPt_
const std::optional< int > maxZ0_
const std::optional< int > minPt_
const std::optional< int > minScalarSumPt_