File indexing completed on 2024-09-12 04:16:40
0001 #ifndef L1Trigger_Phase2L1GT_L1GTSingleCollectionCut_h
0002 #define L1Trigger_Phase2L1GT_L1GTSingleCollectionCut_h
0003
0004 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Utilities/interface/Exception.h"
0007 #include "DataFormats/L1Trigger/interface/P2GTCandidate.h"
0008
0009 #include "L1Trigger/Phase2L1GT/interface/L1GTScales.h"
0010
0011 #include "L1GTOptionalParam.h"
0012
0013 #include <functional>
0014 #include <optional>
0015 #include <cstdint>
0016
0017 namespace l1t {
0018
0019 template <typename T, typename K>
0020 inline std::vector<T> getParamVector(const std::string& name,
0021 const edm::ParameterSet& config,
0022 std::function<T(K)> conv) {
0023 const std::vector<K>& values = config.getParameter<std::vector<K>>(name);
0024 std::vector<T> convertedValues(values.size());
0025 for (std::size_t i = 0; i < values.size(); i++) {
0026 convertedValues[i] = conv(values[i]);
0027 }
0028 return convertedValues;
0029 }
0030
0031 class L1GTSingleCollectionCut {
0032 public:
0033 L1GTSingleCollectionCut(const edm::ParameterSet& config,
0034 const edm::ParameterSet& lutConfig,
0035 const L1GTScales& scales)
0036 : scales_(scales),
0037 tag_(config.getParameter<edm::InputTag>("tag")),
0038 minPt_(getOptionalParam<int, double>(
0039 "minPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
0040 maxPt_(getOptionalParam<int, double>(
0041 "maxPt", config, [&scales](double value) { return scales.to_hw_pT_ceil(value); })),
0042 minEta_(getOptionalParam<int, double>(
0043 "minEta", config, [&scales](double value) { return scales.to_hw_eta_floor(value); })),
0044 maxEta_(getOptionalParam<int, double>(
0045 "maxEta", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
0046 minPhi_(getOptionalParam<int, double>(
0047 "minPhi", config, [&scales](double value) { return scales.to_hw_phi_floor(value); })),
0048 maxPhi_(getOptionalParam<int, double>(
0049 "maxPhi", config, [&scales](double value) { return scales.to_hw_phi_ceil(value); })),
0050 minZ0_(getOptionalParam<int, double>(
0051 "minZ0", config, [&scales](double value) { return scales.to_hw_z0_floor(value); })),
0052 maxZ0_(getOptionalParam<int, double>(
0053 "maxZ0", config, [&scales](double value) { return scales.to_hw_z0_ceil(value); })),
0054 minScalarSumPt_(getOptionalParam<int, double>(
0055 "minScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
0056 maxScalarSumPt_(getOptionalParam<int, double>(
0057 "maxScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT_ceil(value); })),
0058 minQualityScore_(getOptionalParam<unsigned int>("minQualityScore", config)),
0059 maxQualityScore_(getOptionalParam<unsigned int>("maxQualityScore", config)),
0060 qualityFlags_(getOptionalParam<unsigned int>("qualityFlags", config)),
0061 minAbsEta_(getOptionalParam<int, double>(
0062 "minAbsEta", config, [&scales](double value) { return scales.to_hw_eta_floor(value); })),
0063 maxAbsEta_(getOptionalParam<int, double>(
0064 "maxAbsEta", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
0065 minIsolationPt_(getOptionalParam<int, double>(
0066 "minRelIsolationPt", config, [&scales](double value) { return scales.to_hw_isolationPT_floor(value); })),
0067 maxIsolationPt_(getOptionalParam<int, double>(
0068 "maxRelIsolationPt", config, [&scales](double value) { return scales.to_hw_isolationPT_ceil(value); })),
0069 minRelIsolationPt_(getOptionalParam<int, double>(
0070 "minRelIsolationPt",
0071 config,
0072 [&scales](double value) { return scales.to_hw_relative_isolationPT_floor(value); })),
0073 maxRelIsolationPt_(getOptionalParam<int, double>(
0074 "maxRelIsolationPt",
0075 config,
0076 [&scales](double value) { return scales.to_hw_relative_isolationPT_ceil(value); })),
0077 regionsAbsEtaLowerBounds_(getParamVector<int, double>(
0078 "regionsAbsEtaLowerBounds", config, [&scales](double value) { return scales.to_hw_eta_ceil(value); })),
0079 regionsMinPt_(getParamVector<int, double>(
0080 "regionsMinPt", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })),
0081 regionsMaxRelIsolationPt_(getParamVector<int, double>(
0082 "regionsMaxRelIsolationPt",
0083 config,
0084 [&scales](double value) { return scales.to_hw_relative_isolationPT_ceil(value); })),
0085 regionsQualityFlags_(config.getParameter<std::vector<unsigned int>>("regionsQualityFlags")),
0086 minPrimVertDz_(getOptionalParam<int, double>(
0087 "minPrimVertDz", config, [&scales](double value) { return scales.to_hw_z0_floor(value); })),
0088 maxPrimVertDz_(getOptionalParam<int, double>(
0089 "maxPrimVertDz", config, [&scales](double value) { return scales.to_hw_z0_ceil(value); })),
0090 primVertex_(getOptionalParam<unsigned int>("primVertex", config)),
0091 minPtMultiplicityN_(config.getParameter<unsigned int>("minPtMultiplicityN")),
0092 minPtMultiplicityCut_(getOptionalParam<int, double>(
0093 "minPtMultiplicityCut", config, [&scales](double value) { return scales.to_hw_pT_floor(value); })) {}
0094
0095 bool checkObject(const P2GTCandidate& obj) const {
0096 bool result = true;
0097
0098 result &= minPt_ ? (obj.hwPT() > minPt_) : true;
0099 result &= maxPt_ ? (obj.hwPT() < maxPt_) : true;
0100
0101 result &= minEta_ ? (obj.hwEta() > minEta_) : true;
0102 result &= maxEta_ ? (obj.hwEta() < maxEta_) : true;
0103
0104 result &= minPhi_ ? (obj.hwPhi() > minPhi_) : true;
0105 result &= maxPhi_ ? (obj.hwPhi() < maxPhi_) : true;
0106
0107 result &= minZ0_ ? (obj.hwZ0() > minZ0_) : true;
0108 result &= maxZ0_ ? (obj.hwZ0() < maxZ0_) : true;
0109
0110 result &= minAbsEta_ ? (abs(obj.hwEta()) > minAbsEta_) : true;
0111 result &= maxAbsEta_ ? (abs(obj.hwEta()) < maxAbsEta_) : true;
0112
0113 result &= minScalarSumPt_ ? (obj.hwScalarSumPT() > minScalarSumPt_) : true;
0114 result &= maxScalarSumPt_ ? (obj.hwScalarSumPT() < maxScalarSumPt_) : true;
0115
0116 result &= minQualityScore_ ? (obj.hwQualityScore() > minQualityScore_) : true;
0117 result &= maxQualityScore_ ? (obj.hwQualityScore() < maxQualityScore_) : true;
0118 result &= qualityFlags_ ? (obj.hwQualityFlags().to_uint() & qualityFlags_.value()) == qualityFlags_ : true;
0119
0120 result &= minIsolationPt_ ? (obj.hwIsolationPT() > minIsolationPt_) : true;
0121 result &= maxIsolationPt_ ? (obj.hwIsolationPT() < maxIsolationPt_) : true;
0122
0123 result &= minRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION >
0124 static_cast<int64_t>(minRelIsolationPt_.value()) * obj.hwPT().to_int64()
0125 : true;
0126 result &= maxRelIsolationPt_ ? obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
0127 static_cast<int64_t>(maxRelIsolationPt_.value()) * obj.hwPT().to_int64()
0128 : true;
0129
0130 result &= regionsAbsEtaLowerBounds_.empty() ? true : checkEtaDependentCuts(obj);
0131 return result;
0132 }
0133
0134 bool checkCollection(const P2GTCandidateCollection& col) const {
0135 if (minPtMultiplicityN_ > 0 && minPtMultiplicityCut_) {
0136 unsigned int minPtMultiplicity = 0;
0137
0138 for (const P2GTCandidate& obj : col) {
0139 minPtMultiplicity = obj.hwPT() > minPtMultiplicityCut_ ? minPtMultiplicity + 1 : minPtMultiplicity;
0140 }
0141
0142 return minPtMultiplicity >= minPtMultiplicityN_;
0143 }
0144
0145 return true;
0146 }
0147
0148 bool checkPrimaryVertices(const P2GTCandidate& obj, const P2GTCandidateCollection& primVertCol) const {
0149 if (!minPrimVertDz_ && !maxPrimVertDz_) {
0150 return true;
0151 } else {
0152 if (!primVertex_) {
0153 throw cms::Exception("Configuration")
0154 << "When a min/max primary vertex cut is set the primary vertex to cut\n"
0155 << "on has to be specified with with config parameter \'primVertex\'. ";
0156 }
0157 if (primVertex_ < primVertCol.size()) {
0158 uint32_t dZ = abs(obj.hwZ0() - primVertCol[primVertex_.value()].hwZ0());
0159
0160 bool result = true;
0161 result &= minPrimVertDz_ ? dZ > minPrimVertDz_ : true;
0162 result &= maxPrimVertDz_ ? dZ < maxPrimVertDz_ : true;
0163 return result;
0164 } else {
0165 return false;
0166 }
0167 }
0168 }
0169
0170 static void fillPSetDescription(edm::ParameterSetDescription& desc) {
0171 desc.add<edm::InputTag>("tag");
0172 desc.addOptional<double>("minPt");
0173 desc.addOptional<double>("maxPt");
0174 desc.addOptional<double>("minEta");
0175 desc.addOptional<double>("maxEta");
0176 desc.addOptional<double>("minPhi");
0177 desc.addOptional<double>("maxPhi");
0178 desc.addOptional<double>("minZ0");
0179 desc.addOptional<double>("maxZ0");
0180 desc.addOptional<double>("minScalarSumPt");
0181 desc.addOptional<double>("maxScalarSumPt");
0182 desc.addOptional<unsigned int>("minQualityScore");
0183 desc.addOptional<unsigned int>("maxQualityScore");
0184 desc.addOptional<unsigned int>("qualityFlags");
0185 desc.add<std::vector<unsigned int>>("regions", {});
0186 desc.addOptional<double>("minAbsEta");
0187 desc.addOptional<double>("maxAbsEta");
0188 desc.addOptional<double>("minIsolationPt");
0189 desc.addOptional<double>("maxIsolationPt");
0190 desc.addOptional<double>("minRelIsolationPt");
0191 desc.addOptional<double>("maxRelIsolationPt");
0192 desc.add<std::vector<double>>("regionsAbsEtaLowerBounds", {});
0193 desc.add<std::vector<double>>("regionsMinPt", {});
0194 desc.add<std::vector<double>>("regionsMaxRelIsolationPt", {});
0195 desc.add<std::vector<unsigned int>>("regionsQualityFlags", {});
0196 desc.addOptional<double>("minPrimVertDz");
0197 desc.addOptional<double>("maxPrimVertDz");
0198 desc.addOptional<unsigned int>("primVertex");
0199 desc.add<unsigned int>("minPtMultiplicityN", 0);
0200 desc.addOptional<double>("minPtMultiplicityCut");
0201 }
0202
0203 const edm::InputTag& tag() const { return tag_; }
0204
0205 private:
0206 bool checkEtaDependentCuts(const P2GTCandidate& obj) const {
0207 bool res = true;
0208
0209 unsigned int index;
0210 index = atIndex(obj.hwEta());
0211 res &= regionsMinPt_.empty() ? true : obj.hwPT() > regionsMinPt_[index];
0212 res &= regionsMaxRelIsolationPt_.empty()
0213 ? true
0214 : obj.hwIsolationPT().to_int64() << L1GTScales::RELATIVE_ISOLATION_RESOLUTION <
0215 static_cast<int64_t>(regionsMaxRelIsolationPt_[index]) * obj.hwPT().to_int64();
0216 res &= regionsQualityFlags_.empty()
0217 ? true
0218 : (obj.hwQualityFlags().to_uint() & regionsQualityFlags_[index]) == regionsQualityFlags_[index];
0219 return res;
0220 }
0221
0222 unsigned int atIndex(int objeta) const {
0223
0224
0225
0226 for (unsigned int i = 0; i < regionsAbsEtaLowerBounds_.size() - 1; i++) {
0227 if (std::abs(objeta) >= regionsAbsEtaLowerBounds_[i] && std::abs(objeta) < regionsAbsEtaLowerBounds_[i + 1]) {
0228 return i;
0229 }
0230 }
0231 return regionsAbsEtaLowerBounds_.size() - 1;
0232 }
0233
0234 private:
0235 const L1GTScales scales_;
0236 const edm::InputTag tag_;
0237 const std::optional<int> minPt_;
0238 const std::optional<int> maxPt_;
0239 const std::optional<int> minEta_;
0240 const std::optional<int> maxEta_;
0241 const std::optional<int> minPhi_;
0242 const std::optional<int> maxPhi_;
0243 const std::optional<int> minZ0_;
0244 const std::optional<int> maxZ0_;
0245 const std::optional<int> minScalarSumPt_;
0246 const std::optional<int> maxScalarSumPt_;
0247 const std::optional<unsigned int> minQualityScore_;
0248 const std::optional<unsigned int> maxQualityScore_;
0249 const std::optional<unsigned int> qualityFlags_;
0250 const std::optional<int> minAbsEta_;
0251 const std::optional<int> maxAbsEta_;
0252 const std::optional<int> minIsolationPt_;
0253 const std::optional<int> maxIsolationPt_;
0254 const std::optional<int> minRelIsolationPt_;
0255 const std::optional<int> maxRelIsolationPt_;
0256 const std::vector<int> regionsAbsEtaLowerBounds_;
0257 const std::vector<int> regionsMinPt_;
0258 const std::vector<int> regionsMaxRelIsolationPt_;
0259 const std::vector<unsigned int> regionsQualityFlags_;
0260 const std::optional<int> minPrimVertDz_;
0261 const std::optional<int> maxPrimVertDz_;
0262 const std::optional<unsigned int> primVertex_;
0263 const unsigned int minPtMultiplicityN_;
0264 const std::optional<int> minPtMultiplicityCut_;
0265 };
0266
0267 }
0268
0269 #endif