File indexing completed on 2023-10-25 09:55:51
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 "DataFormats/L1Trigger/interface/P2GTCandidate.h"
0007
0008 #include "L1Trigger/Phase2L1GT/interface/L1GTScales.h"
0009
0010 #include "L1GTOptionalParam.h"
0011
0012 #include <functional>
0013 #include <optional>
0014
0015 namespace l1t {
0016
0017 template <typename T, typename K>
0018 inline std::vector<T> getParamVector(const std::string& name,
0019 const edm::ParameterSet& config,
0020 std::function<T(K)> conv) {
0021 const std::vector<K>& values = config.getParameter<std::vector<K>>(name);
0022 std::vector<T> convertedValues(values.size());
0023 for (std::size_t i = 0; i < values.size(); i++) {
0024 convertedValues[i] = conv(values[i]);
0025 }
0026 return convertedValues;
0027 }
0028
0029 class L1GTSingleCollectionCut {
0030 public:
0031 L1GTSingleCollectionCut(const edm::ParameterSet& config,
0032 const edm::ParameterSet& lutConfig,
0033 const L1GTScales& scales)
0034 : scales_(scales),
0035 tag_(config.getParameter<edm::InputTag>("tag")),
0036 minPt_(getOptionalParam<int, double>(
0037 "minPt", config, [&scales](double value) { return scales.to_hw_pT(value); })),
0038 maxPt_(getOptionalParam<int, double>(
0039 "maxPt", config, [&scales](double value) { return scales.to_hw_pT(value); })),
0040 minEta_(getOptionalParam<int, double>(
0041 "minEta", config, [&scales](double value) { return scales.to_hw_eta(value); })),
0042 maxEta_(getOptionalParam<int, double>(
0043 "maxEta", config, [&scales](double value) { return scales.to_hw_eta(value); })),
0044 minPhi_(getOptionalParam<int, double>(
0045 "minPhi", config, [&scales](double value) { return scales.to_hw_phi(value); })),
0046 maxPhi_(getOptionalParam<int, double>(
0047 "maxPhi", config, [&scales](double value) { return scales.to_hw_phi(value); })),
0048 minZ0_(getOptionalParam<int, double>(
0049 "minZ0", config, [&scales](double value) { return scales.to_hw_z0(value); })),
0050 maxZ0_(getOptionalParam<int, double>(
0051 "maxZ0", config, [&scales](double value) { return scales.to_hw_z0(value); })),
0052 minScalarSumPt_(getOptionalParam<int, double>(
0053 "minScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT(value); })),
0054 maxScalarSumPt_(getOptionalParam<int, double>(
0055 "maxScalarSumPt", config, [&scales](double value) { return scales.to_hw_pT(value); })),
0056 qual_(config.getParameter<std::vector<unsigned int>>("qual")),
0057 minAbsEta_(getOptionalParam<int, double>(
0058 "minAbsEta", config, [&scales](double value) { return scales.to_hw_eta(value); })),
0059 maxAbsEta_(getOptionalParam<int, double>(
0060 "maxAbsEta", config, [&scales](double value) { return scales.to_hw_eta(value); })),
0061 maxIso_(getOptionalParam<int, double>(
0062 "maxIso", config, [&scales](double value) { return scales.to_hw_isolation(value); })),
0063 minHwIso_(getOptionalParam<int>("minHwIso", config)),
0064 regionsAbsEtaLowerBounds_(getParamVector<int, double>(
0065 "regionsAbsEtaLowerBounds", config, [&scales](double value) { return scales.to_hw_eta(value); })),
0066 regionsMinPt_(getParamVector<int, double>(
0067 "regionsMinPt", config, [&scales](double value) { return scales.to_hw_pT(value); })),
0068 regionsMaxIso_(getParamVector<int, double>(
0069 "regionsMaxIso", config, [&scales](double value) { return scales.to_hw_isolation(value); })),
0070 regionsQual_(config.getParameter<std::vector<unsigned int>>("regionsQual")) {}
0071
0072 bool checkEtadependentcuts(const P2GTCandidate& obj) const {
0073 bool res = true;
0074
0075 unsigned int index;
0076 index = atIndex(obj.hwEta());
0077 res &= regionsMinPt_.empty() ? true : obj.hwPT() > regionsMinPt_[index];
0078 res &= regionsMaxIso_.empty()
0079 ? true
0080 : obj.hwIso().to_int() << scales_.isolation_shift() < regionsMaxIso_[index] * obj.hwPT().to_int();
0081 res &= regionsQual_.empty() ? true : (obj.hwQual().to_uint() & regionsQual_[index]) == regionsQual_[index];
0082 return res;
0083 }
0084
0085 unsigned int atIndex(int objeta) const {
0086
0087
0088
0089 for (unsigned int i = 0; i < regionsAbsEtaLowerBounds_.size() - 1; i++) {
0090 if (std::abs(objeta) >= regionsAbsEtaLowerBounds_[i] && std::abs(objeta) < regionsAbsEtaLowerBounds_[i + 1]) {
0091 return i;
0092 }
0093 }
0094 return regionsAbsEtaLowerBounds_.size() - 1;
0095 }
0096
0097 bool checkObject(const P2GTCandidate& obj) const {
0098 bool result = true;
0099
0100 result &= minPt_ ? (obj.hwPT() > minPt_) : true;
0101 result &= maxPt_ ? (obj.hwPT() < maxPt_) : true;
0102
0103 result &= minEta_ ? (obj.hwEta() > minEta_) : true;
0104 result &= maxEta_ ? (obj.hwEta() < maxEta_) : true;
0105
0106 result &= minPhi_ ? (obj.hwPhi() > minPhi_) : true;
0107 result &= maxPhi_ ? (obj.hwPhi() < maxPhi_) : true;
0108
0109 result &= minZ0_ ? (obj.hwZ0() > minZ0_) : true;
0110 result &= maxZ0_ ? (obj.hwZ0() < maxZ0_) : true;
0111
0112 result &= minAbsEta_ ? (abs(obj.hwEta()) > minAbsEta_) : true;
0113 result &= maxAbsEta_ ? (abs(obj.hwEta()) < maxAbsEta_) : true;
0114
0115 result &= minScalarSumPt_ ? (obj.hwSca_sum() > minScalarSumPt_) : true;
0116 result &= maxScalarSumPt_ ? (obj.hwSca_sum() < minScalarSumPt_) : true;
0117
0118 result &= qual_.empty() ? true : std::find(qual_.begin(), qual_.end(), obj.hwQual().to_uint()) != qual_.end();
0119 result &=
0120 maxIso_ ? obj.hwIso().to_int() << scales_.isolation_shift() < maxIso_.value() * obj.hwPT().to_int() : true;
0121
0122 result &= minHwIso_ ? (obj.hwIso() > minHwIso_) : true;
0123 result &= regionsAbsEtaLowerBounds_.empty() ? true : checkEtadependentcuts(obj);
0124 return result;
0125 }
0126
0127 static void fillPSetDescription(edm::ParameterSetDescription& desc) {
0128 desc.add<edm::InputTag>("tag");
0129 desc.addOptional<double>("minPt");
0130 desc.addOptional<double>("maxPt");
0131 desc.addOptional<double>("minEta");
0132 desc.addOptional<double>("maxEta");
0133 desc.addOptional<double>("minPhi");
0134 desc.addOptional<double>("maxPhi");
0135 desc.addOptional<double>("minZ0");
0136 desc.addOptional<double>("maxZ0");
0137 desc.addOptional<double>("minScalarSumPt");
0138 desc.addOptional<double>("maxScalarSumPt");
0139 desc.add<std::vector<unsigned int>>("qual", {});
0140 desc.addOptional<double>("minAbsEta");
0141 desc.addOptional<double>("maxAbsEta");
0142 desc.addOptional<double>("maxIso");
0143 desc.addOptional<int>("minHwIso");
0144 desc.add<std::vector<double>>("regionsAbsEtaLowerBounds", {});
0145 desc.add<std::vector<double>>("regionsMinPt", {});
0146 desc.add<std::vector<double>>("regionsMaxIso", {});
0147 desc.add<std::vector<unsigned int>>("regionsQual", {});
0148 }
0149
0150 const edm::InputTag& tag() const { return tag_; }
0151
0152 private:
0153 const L1GTScales scales_;
0154 const edm::InputTag tag_;
0155 const std::optional<int> minPt_;
0156 const std::optional<int> maxPt_;
0157 const std::optional<int> minEta_;
0158 const std::optional<int> maxEta_;
0159 const std::optional<int> minPhi_;
0160 const std::optional<int> maxPhi_;
0161 const std::optional<int> minZ0_;
0162 const std::optional<int> maxZ0_;
0163 const std::optional<int> minScalarSumPt_;
0164 const std::optional<int> maxScalarSumPt_;
0165 const std::vector<unsigned int> qual_;
0166 const std::optional<int> minAbsEta_;
0167 const std::optional<int> maxAbsEta_;
0168 const std::optional<int> maxIso_;
0169 const std::optional<int> minHwIso_;
0170 const std::vector<int> regionsAbsEtaLowerBounds_;
0171 const std::vector<int> regionsMinPt_;
0172 const std::vector<int> regionsMaxIso_;
0173 const std::vector<unsigned int> regionsQual_;
0174 };
0175
0176 }
0177
0178 #endif