CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
l1t::L1GTSingleCollectionCut Class Reference

#include <L1GTSingleCollectionCut.h>

Public Member Functions

bool checkCollection (const P2GTCandidateCollection &col) const
 
bool checkObject (const P2GTCandidate &obj) const
 
bool checkPrimaryVertices (const P2GTCandidate &obj, const P2GTCandidateCollection &primVertCol) const
 
 L1GTSingleCollectionCut (const edm::ParameterSet &config, const edm::ParameterSet &lutConfig, const L1GTScales &scales)
 
const edm::InputTagtag () const
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &desc)
 

Private Member Functions

unsigned int atIndex (int objeta) const
 
bool checkEtaDependentCuts (const P2GTCandidate &obj) const
 

Private Attributes

const std::optional< int > maxAbsEta_
 
const std::optional< int > maxEta_
 
const std::optional< int > maxIsolationPt_
 
const std::optional< int > maxPhi_
 
const std::optional< int > maxPrimVertDz_
 
const std::optional< int > maxPt_
 
const std::optional< unsigned int > maxQualityScore_
 
const std::optional< int > maxRelIsolationPt_
 
const std::optional< int > maxScalarSumPt_
 
const std::optional< int > maxZ0_
 
const std::optional< int > minAbsEta_
 
const std::optional< int > minEta_
 
const std::optional< int > minIsolationPt_
 
const std::optional< int > minPhi_
 
const std::optional< int > minPrimVertDz_
 
const std::optional< int > minPt_
 
const std::optional< int > minPtMultiplicityCut_
 
const unsigned int minPtMultiplicityN_
 
const std::optional< unsigned int > minQualityScore_
 
const std::optional< int > minRelIsolationPt_
 
const std::optional< int > minScalarSumPt_
 
const std::optional< int > minZ0_
 
const std::optional< unsigned int > primVertex_
 
const std::optional< unsigned int > qualityFlags_
 
const std::vector< int > regionsAbsEtaLowerBounds_
 
const std::vector< int > regionsMaxRelIsolationPt_
 
const std::vector< int > regionsMinPt_
 
const std::vector< unsigned int > regionsQualityFlags_
 
const L1GTScales scales_
 
const edm::InputTag tag_
 

Detailed Description

Definition at line 31 of file L1GTSingleCollectionCut.h.

Constructor & Destructor Documentation

◆ L1GTSingleCollectionCut()

l1t::L1GTSingleCollectionCut::L1GTSingleCollectionCut ( const edm::ParameterSet config,
const edm::ParameterSet lutConfig,
const L1GTScales scales 
)
inline

Definition at line 33 of file L1GTSingleCollectionCut.h.

References l1tGTDoubleObjectCond_cfi::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  minPtMultiplicityN_(config.getParameter<unsigned int>("minPtMultiplicityN")),
92  minPtMultiplicityCut_(getOptionalParam<int, double>(
93  "minPtMultiplicityCut", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })) {}
const std::optional< int > maxAbsEta_
const std::optional< unsigned int > maxQualityScore_
const std::optional< int > maxRelIsolationPt_
const std::vector< int > regionsMaxRelIsolationPt_
const std::optional< int > minPhi_
const std::optional< unsigned int > primVertex_
const std::optional< int > minPrimVertDz_
const std::vector< unsigned int > regionsQualityFlags_
const std::optional< int > maxEta_
const std::optional< int > maxPhi_
const std::vector< int > regionsAbsEtaLowerBounds_
const std::optional< int > maxScalarSumPt_
Definition: config.py:1
const std::optional< unsigned int > qualityFlags_
const std::optional< int > minAbsEta_
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_
const std::optional< int > minPtMultiplicityCut_
Definition: value.py:1
const std::optional< unsigned int > minQualityScore_
const std::vector< int > regionsMinPt_
const std::optional< int > minIsolationPt_
const std::optional< int > maxZ0_
const std::optional< int > minPt_
const std::optional< int > minScalarSumPt_

Member Function Documentation

◆ atIndex()

unsigned int l1t::L1GTSingleCollectionCut::atIndex ( int  objeta) const
inlineprivate

Definition at line 222 of file L1GTSingleCollectionCut.h.

References funct::abs(), mps_fire::i, and regionsAbsEtaLowerBounds_.

Referenced by checkEtaDependentCuts().

222  {
223  // Function that checks at which index the eta of the object is
224  // If object abs(eta) < regionsAbsEtaLowerBounds_[0] the function returns the last index,
225  // Might be undesired behaviour
226  for (unsigned int i = 0; i < regionsAbsEtaLowerBounds_.size() - 1; i++) {
227  if (std::abs(objeta) >= regionsAbsEtaLowerBounds_[i] && std::abs(objeta) < regionsAbsEtaLowerBounds_[i + 1]) {
228  return i;
229  }
230  }
231  return regionsAbsEtaLowerBounds_.size() - 1;
232  }
const std::vector< int > regionsAbsEtaLowerBounds_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ checkCollection()

bool l1t::L1GTSingleCollectionCut::checkCollection ( const P2GTCandidateCollection col) const
inline

Definition at line 134 of file L1GTSingleCollectionCut.h.

References cuy::col, minPtMultiplicityCut_, minPtMultiplicityN_, and getGTfromDQMFile::obj.

Referenced by L1GTSingleObjectCond::filter(), L1GTTripleObjectCond::filter(), L1GTQuadObjectCond::filter(), and L1GTDoubleObjectCond::filter().

134  {
136  unsigned int minPtMultiplicity = 0;
137 
138  for (const P2GTCandidate& obj : col) {
139  minPtMultiplicity = obj.hwPT() > minPtMultiplicityCut_ ? minPtMultiplicity + 1 : minPtMultiplicity;
140  }
141 
142  return minPtMultiplicity >= minPtMultiplicityN_;
143  }
144 
145  return true;
146  }
const std::optional< int > minPtMultiplicityCut_
col
Definition: cuy.py:1009

◆ checkEtaDependentCuts()

bool l1t::L1GTSingleCollectionCut::checkEtaDependentCuts ( const P2GTCandidate obj) const
inlineprivate

Definition at line 206 of file L1GTSingleCollectionCut.h.

References atIndex(), getGTfromDQMFile::obj, regionsMaxRelIsolationPt_, regionsMinPt_, regionsQualityFlags_, l1t::L1GTScales::RELATIVE_ISOLATION_RESOLUTION, and funct::true.

Referenced by checkObject().

206  {
207  bool res = true;
208 
209  unsigned int index;
210  index = atIndex(obj.hwEta());
211  res &= regionsMinPt_.empty() ? true : obj.hwPT() > regionsMinPt_[index];
212  res &= regionsMaxRelIsolationPt_.empty()
213  ? true
214  : obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
215  static_cast<int64_t>(regionsMaxRelIsolationPt_[index]) * obj.hwPT().to_int64();
216  res &= regionsQualityFlags_.empty()
217  ? true
218  : (obj.hwQualityFlags().to_uint() & regionsQualityFlags_[index]) == regionsQualityFlags_[index];
219  return res;
220  }
const std::vector< int > regionsMaxRelIsolationPt_
const std::vector< unsigned int > regionsQualityFlags_
Definition: Electron.h:6
unsigned int atIndex(int objeta) const
static constexpr int RELATIVE_ISOLATION_RESOLUTION
Definition: L1GTScales.h:12
const std::vector< int > regionsMinPt_

◆ checkObject()

bool l1t::L1GTSingleCollectionCut::checkObject ( const P2GTCandidate obj) const
inline

Definition at line 95 of file L1GTSingleCollectionCut.h.

References funct::abs(), checkEtaDependentCuts(), maxAbsEta_, maxEta_, maxIsolationPt_, maxPhi_, maxPt_, maxQualityScore_, maxRelIsolationPt_, maxScalarSumPt_, maxZ0_, minAbsEta_, minEta_, minIsolationPt_, minPhi_, minPt_, minQualityScore_, minRelIsolationPt_, minScalarSumPt_, minZ0_, getGTfromDQMFile::obj, qualityFlags_, regionsAbsEtaLowerBounds_, l1t::L1GTScales::RELATIVE_ISOLATION_RESOLUTION, mps_fire::result, and funct::true.

Referenced by L1GTSingleObjectCond::filter(), L1GTTripleObjectCond::filter(), L1GTQuadObjectCond::filter(), and L1GTDoubleObjectCond::filter().

95  {
96  bool result = true;
97 
98  result &= minPt_ ? (obj.hwPT() > minPt_) : true;
99  result &= maxPt_ ? (obj.hwPT() < maxPt_) : true;
100 
101  result &= minEta_ ? (obj.hwEta() > minEta_) : true;
102  result &= maxEta_ ? (obj.hwEta() < maxEta_) : true;
103 
104  result &= minPhi_ ? (obj.hwPhi() > minPhi_) : true;
105  result &= maxPhi_ ? (obj.hwPhi() < maxPhi_) : true;
106 
107  result &= minZ0_ ? (obj.hwZ0() > minZ0_) : true;
108  result &= maxZ0_ ? (obj.hwZ0() < maxZ0_) : true;
109 
110  result &= minAbsEta_ ? (abs(obj.hwEta()) > minAbsEta_) : true;
111  result &= maxAbsEta_ ? (abs(obj.hwEta()) < maxAbsEta_) : true;
112 
113  result &= minScalarSumPt_ ? (obj.hwScalarSumPT() > minScalarSumPt_) : true;
114  result &= maxScalarSumPt_ ? (obj.hwScalarSumPT() < maxScalarSumPt_) : true;
115 
116  result &= minQualityScore_ ? (obj.hwQualityScore() > minQualityScore_) : true;
117  result &= maxQualityScore_ ? (obj.hwQualityScore() < maxQualityScore_) : true;
118  result &= qualityFlags_ ? (obj.hwQualityFlags().to_uint() & qualityFlags_.value()) == qualityFlags_ : true;
119 
120  result &= minIsolationPt_ ? (obj.hwIsolationPT() > minIsolationPt_) : true;
121  result &= maxIsolationPt_ ? (obj.hwIsolationPT() < maxIsolationPt_) : true;
122 
123  result &= minRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION >
124  static_cast<int64_t>(minRelIsolationPt_.value()) * obj.hwPT().to_int64()
125  : true;
126  result &= maxRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
127  static_cast<int64_t>(maxRelIsolationPt_.value()) * obj.hwPT().to_int64()
128  : true;
129 
131  return result;
132  }
const std::optional< int > maxAbsEta_
const std::optional< unsigned int > maxQualityScore_
const std::optional< int > maxRelIsolationPt_
bool checkEtaDependentCuts(const P2GTCandidate &obj) const
const std::optional< int > minPhi_
const std::optional< int > maxEta_
const std::optional< int > maxPhi_
const std::vector< int > regionsAbsEtaLowerBounds_
const std::optional< int > maxScalarSumPt_
const std::optional< unsigned int > qualityFlags_
const std::optional< int > minAbsEta_
const std::optional< int > minRelIsolationPt_
const std::optional< int > minEta_
const std::optional< int > maxPt_
const std::optional< int > minZ0_
const std::optional< int > maxIsolationPt_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const std::optional< unsigned int > minQualityScore_
static constexpr int RELATIVE_ISOLATION_RESOLUTION
Definition: L1GTScales.h:12
const std::optional< int > minIsolationPt_
const std::optional< int > maxZ0_
const std::optional< int > minPt_
const std::optional< int > minScalarSumPt_

◆ checkPrimaryVertices()

bool l1t::L1GTSingleCollectionCut::checkPrimaryVertices ( const P2GTCandidate obj,
const P2GTCandidateCollection primVertCol 
) const
inline

Definition at line 148 of file L1GTSingleCollectionCut.h.

References funct::abs(), Exception, maxPrimVertDz_, minPrimVertDz_, getGTfromDQMFile::obj, primVertex_, and mps_fire::result.

Referenced by L1GTSingleObjectCond::filter(), L1GTTripleObjectCond::filter(), L1GTQuadObjectCond::filter(), and L1GTDoubleObjectCond::filter().

148  {
149  if (!minPrimVertDz_ && !maxPrimVertDz_) {
150  return true;
151  } else {
152  if (!primVertex_) {
153  throw cms::Exception("Configuration")
154  << "When a min/max primary vertex cut is set the primary vertex to cut\n"
155  << "on has to be specified with with config parameter \'primVertex\'. ";
156  }
157  if (primVertex_ < primVertCol.size()) {
158  uint32_t dZ = abs(obj.hwZ0() - primVertCol[primVertex_.value()].hwZ0());
159 
160  bool result = true;
161  result &= minPrimVertDz_ ? dZ > minPrimVertDz_ : true;
162  result &= maxPrimVertDz_ ? dZ < maxPrimVertDz_ : true;
163  return result;
164  } else {
165  return false;
166  }
167  }
168  }
const std::optional< unsigned int > primVertex_
const std::optional< int > minPrimVertDz_
const std::optional< int > maxPrimVertDz_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ fillPSetDescription()

static void l1t::L1GTSingleCollectionCut::fillPSetDescription ( edm::ParameterSetDescription desc)
inlinestatic

Definition at line 170 of file L1GTSingleCollectionCut.h.

References submitPVResolutionJobs::desc.

170  {
171  desc.add<edm::InputTag>("tag");
172  desc.addOptional<double>("minPt");
173  desc.addOptional<double>("maxPt");
174  desc.addOptional<double>("minEta");
175  desc.addOptional<double>("maxEta");
176  desc.addOptional<double>("minPhi");
177  desc.addOptional<double>("maxPhi");
178  desc.addOptional<double>("minZ0");
179  desc.addOptional<double>("maxZ0");
180  desc.addOptional<double>("minScalarSumPt");
181  desc.addOptional<double>("maxScalarSumPt");
182  desc.addOptional<unsigned int>("minQualityScore");
183  desc.addOptional<unsigned int>("maxQualityScore");
184  desc.addOptional<unsigned int>("qualityFlags");
185  desc.add<std::vector<unsigned int>>("regions", {});
186  desc.addOptional<double>("minAbsEta");
187  desc.addOptional<double>("maxAbsEta");
188  desc.addOptional<double>("minIsolationPt");
189  desc.addOptional<double>("maxIsolationPt");
190  desc.addOptional<double>("minRelIsolationPt");
191  desc.addOptional<double>("maxRelIsolationPt");
192  desc.add<std::vector<double>>("regionsAbsEtaLowerBounds", {});
193  desc.add<std::vector<double>>("regionsMinPt", {});
194  desc.add<std::vector<double>>("regionsMaxRelIsolationPt", {});
195  desc.add<std::vector<unsigned int>>("regionsQualityFlags", {});
196  desc.addOptional<double>("minPrimVertDz");
197  desc.addOptional<double>("maxPrimVertDz");
198  desc.addOptional<unsigned int>("primVertex");
199  desc.add<unsigned int>("minPtMultiplicityN", 0);
200  desc.addOptional<double>("minPtMultiplicityCut");
201  }

◆ tag()

const edm::InputTag& l1t::L1GTSingleCollectionCut::tag ( ) const
inline

Member Data Documentation

◆ maxAbsEta_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxAbsEta_
private

Definition at line 251 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxEta_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxEta_
private

Definition at line 240 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxIsolationPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxIsolationPt_
private

Definition at line 253 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxPhi_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxPhi_
private

Definition at line 242 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxPrimVertDz_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxPrimVertDz_
private

Definition at line 261 of file L1GTSingleCollectionCut.h.

Referenced by checkPrimaryVertices().

◆ maxPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxPt_
private

Definition at line 238 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxQualityScore_

const std::optional<unsigned int> l1t::L1GTSingleCollectionCut::maxQualityScore_
private

Definition at line 248 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxRelIsolationPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxRelIsolationPt_
private

Definition at line 255 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxScalarSumPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxScalarSumPt_
private

Definition at line 246 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ maxZ0_

const std::optional<int> l1t::L1GTSingleCollectionCut::maxZ0_
private

Definition at line 244 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minAbsEta_

const std::optional<int> l1t::L1GTSingleCollectionCut::minAbsEta_
private

Definition at line 250 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minEta_

const std::optional<int> l1t::L1GTSingleCollectionCut::minEta_
private

Definition at line 239 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minIsolationPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::minIsolationPt_
private

Definition at line 252 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minPhi_

const std::optional<int> l1t::L1GTSingleCollectionCut::minPhi_
private

Definition at line 241 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minPrimVertDz_

const std::optional<int> l1t::L1GTSingleCollectionCut::minPrimVertDz_
private

Definition at line 260 of file L1GTSingleCollectionCut.h.

Referenced by checkPrimaryVertices().

◆ minPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::minPt_
private

Definition at line 237 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minPtMultiplicityCut_

const std::optional<int> l1t::L1GTSingleCollectionCut::minPtMultiplicityCut_
private

Definition at line 264 of file L1GTSingleCollectionCut.h.

Referenced by checkCollection().

◆ minPtMultiplicityN_

const unsigned int l1t::L1GTSingleCollectionCut::minPtMultiplicityN_
private

Definition at line 263 of file L1GTSingleCollectionCut.h.

Referenced by checkCollection().

◆ minQualityScore_

const std::optional<unsigned int> l1t::L1GTSingleCollectionCut::minQualityScore_
private

Definition at line 247 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minRelIsolationPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::minRelIsolationPt_
private

Definition at line 254 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minScalarSumPt_

const std::optional<int> l1t::L1GTSingleCollectionCut::minScalarSumPt_
private

Definition at line 245 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ minZ0_

const std::optional<int> l1t::L1GTSingleCollectionCut::minZ0_
private

Definition at line 243 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ primVertex_

const std::optional<unsigned int> l1t::L1GTSingleCollectionCut::primVertex_
private

Definition at line 262 of file L1GTSingleCollectionCut.h.

Referenced by checkPrimaryVertices().

◆ qualityFlags_

const std::optional<unsigned int> l1t::L1GTSingleCollectionCut::qualityFlags_
private

Definition at line 249 of file L1GTSingleCollectionCut.h.

Referenced by checkObject().

◆ regionsAbsEtaLowerBounds_

const std::vector<int> l1t::L1GTSingleCollectionCut::regionsAbsEtaLowerBounds_
private

Definition at line 256 of file L1GTSingleCollectionCut.h.

Referenced by atIndex(), and checkObject().

◆ regionsMaxRelIsolationPt_

const std::vector<int> l1t::L1GTSingleCollectionCut::regionsMaxRelIsolationPt_
private

Definition at line 258 of file L1GTSingleCollectionCut.h.

Referenced by checkEtaDependentCuts().

◆ regionsMinPt_

const std::vector<int> l1t::L1GTSingleCollectionCut::regionsMinPt_
private

Definition at line 257 of file L1GTSingleCollectionCut.h.

Referenced by checkEtaDependentCuts().

◆ regionsQualityFlags_

const std::vector<unsigned int> l1t::L1GTSingleCollectionCut::regionsQualityFlags_
private

Definition at line 259 of file L1GTSingleCollectionCut.h.

Referenced by checkEtaDependentCuts().

◆ scales_

const L1GTScales l1t::L1GTSingleCollectionCut::scales_
private

Definition at line 235 of file L1GTSingleCollectionCut.h.

◆ tag_

const edm::InputTag l1t::L1GTSingleCollectionCut::tag_
private

Definition at line 236 of file L1GTSingleCollectionCut.h.

Referenced by tag().