File indexing completed on 2025-02-13 02:58:34
0001 #include "Validation/HLTrigger/interface/HLTGenResHistColl.h"
0002
0003 namespace {
0004
0005 std::pair<trigger::size_type, trigger::size_type> getTrigObjIndicesOfCollection(
0006 const trigger::TriggerEvent& triggerEvent, const std::string& collectionName, const std::string& hltProcessName) {
0007 trigger::size_type begin = 0;
0008 trigger::size_type end = 0;
0009 for (trigger::size_type collNr = 0; collNr < triggerEvent.sizeCollections(); ++collNr) {
0010 std::string collection = triggerEvent.collectionTag(collNr).label();
0011 if (collection == collectionName) {
0012 if (hltProcessName.empty() || triggerEvent.collectionTag(collNr).process() == hltProcessName) {
0013 if (collNr == 0) {
0014 begin = 0;
0015 } else {
0016 begin = triggerEvent.collectionKey(collNr - 1);
0017 }
0018 end = triggerEvent.collectionKey(collNr);
0019 break;
0020 }
0021 }
0022 }
0023 return std::make_pair(begin, end);
0024 }
0025
0026 bool passesFilter(trigger::size_type objKey,
0027 const trigger::TriggerEvent& trigEvent,
0028 const std::string& filterName,
0029 const std::string& hltProcessName) {
0030 edm::InputTag filterTag(filterName, "", hltProcessName);
0031 trigger::size_type filterIndex = trigEvent.filterIndex(filterTag);
0032 if (filterIndex < trigEvent.sizeFilters()) {
0033 const trigger::Keys& trigKeys = trigEvent.filterKeys(filterIndex);
0034 return std::find(trigKeys.begin(), trigKeys.end(), objKey) != trigKeys.end();
0035 } else {
0036 return false;
0037 }
0038 }
0039 template <typename InType, typename OutType>
0040 std::vector<OutType> convertVec(const std::vector<InType>& inVec) {
0041 std::vector<OutType> outVec;
0042 outVec.reserve(inVec.size());
0043 for (const auto& inVal : inVec) {
0044 outVec.push_back(static_cast<OutType>(inVal));
0045 }
0046 return outVec;
0047 }
0048 }
0049
0050
0051 HLTGenResHistColl::HLTGenResHistColl(edm::ParameterSet filterCollConfig, std::string hltProcessName)
0052 : hltProcessName_(hltProcessName) {
0053 objType_ = filterCollConfig.getParameter<std::string>("objType");
0054 isEventLevelVariable_ = objType_ == "MET" || objType_ == "AK4HT" || objType_ == "AK8HT";
0055 filters_ = filterCollConfig.getParameter<std::vector<std::string>>("filterNames");
0056 andFilters_ = filterCollConfig.getParameter<bool>("andFilters");
0057 collectionName_ = filterCollConfig.getParameter<std::string>("collectionName");
0058 dR2limit_ = filterCollConfig.getParameter<double>("dR2limit");
0059 histConfigs_ = filterCollConfig.getParameter<std::vector<edm::ParameterSet>>("histConfigs");
0060 histConfigs2D_ = filterCollConfig.getParameter<std::vector<edm::ParameterSet>>("histConfigs2D");
0061 histNamePrefix_ = "resHist";
0062 separator_ = "__";
0063 }
0064
0065 edm::ParameterSetDescription HLTGenResHistColl::makePSetDescription() {
0066 edm::ParameterSetDescription desc;
0067 desc.add<std::string>("objType", "");
0068 desc.add<std::vector<std::string>>("filterNames", {});
0069 desc.add<bool>("andFilters", true);
0070 desc.add<std::string>("collectionName", "");
0071 desc.add<double>("dR2limit", 0.1);
0072 std::vector<edm::ParameterSet> histConfigDefaults;
0073
0074 edm::ParameterSetDescription histConfig;
0075 histConfig.add<std::string>("vsVar");
0076 histConfig.add<std::string>("resVar");
0077 histConfig.add<std::vector<double>>("vsBinLowEdges");
0078 histConfig.add<std::vector<double>>("resBinLowEdges");
0079 histConfig.addVPSet(
0080 "rangeCuts", VarRangeCut<HLTGenValObject>::makePSetDescription(), std::vector<edm::ParameterSet>());
0081
0082 edm::ParameterSet histConfigDefault0;
0083 histConfigDefault0.addParameter<std::string>("resVar", "ptRes");
0084 histConfigDefault0.addParameter<std::string>("vsVar", "pt");
0085 std::vector<double> defaultPtBinning{0, 5, 10, 12.5, 15, 17.5, 20, 22.5, 25, 30, 35, 40,
0086 45, 50, 60, 80, 100, 150, 200, 250, 300, 350, 400};
0087 histConfigDefault0.addParameter<std::vector<double>>("vsBinLowEdges", defaultPtBinning);
0088 std::vector<double> defaultResBinning;
0089 defaultResBinning.reserve(151);
0090 for (int i = 0; i < 151; i++) {
0091 defaultResBinning.push_back(i * 0.01);
0092 }
0093 histConfigDefault0.addParameter<std::vector<double>>("resBinLowEdges", defaultResBinning);
0094
0095 histConfigDefaults.push_back(histConfigDefault0);
0096 desc.addVPSet("histConfigs", histConfig, histConfigDefaults);
0097
0098
0099 edm::ParameterSetDescription histConfig2D;
0100 histConfig2D.add<std::string>("vsVarX");
0101 histConfig2D.add<std::string>("vsVarY");
0102 histConfig2D.add<std::vector<double>>("binLowEdgesX");
0103 histConfig2D.add<std::vector<double>>("binLowEdgesY");
0104 histConfig2D.addVPSet(
0105 "rangeCuts", VarRangeCut<HLTGenValObject>::makePSetDescription(), std::vector<edm::ParameterSet>());
0106
0107 std::vector<edm::ParameterSet> histConfigDefaults2D;
0108 desc.addVPSet("histConfigs2D", histConfig2D, histConfigDefaults2D);
0109
0110 return desc;
0111 }
0112
0113 bool HLTGenResHistColl::passFilterSelection(trigger::size_type objKey,
0114 const trigger::TriggerEvent& triggerEvent) const {
0115 if (andFilters_) {
0116 for (const auto& filter : filters_) {
0117 if (!passesFilter(objKey, triggerEvent, filter, hltProcessName_)) {
0118 return false;
0119 }
0120 }
0121 return true;
0122 } else {
0123 for (const auto& filter : filters_) {
0124 if (passesFilter(objKey, triggerEvent, filter, hltProcessName_)) {
0125 return true;
0126 }
0127 }
0128 return false;
0129 }
0130 }
0131
0132
0133 void HLTGenResHistColl::bookHists(DQMStore::IBooker& iBooker) {
0134 for (const auto& histConfig : histConfigs_)
0135 book1D(iBooker, histConfig);
0136 for (const auto& histConfig : histConfigs2D_)
0137 book2D(iBooker, histConfig);
0138 }
0139
0140
0141 void HLTGenResHistColl::fillHists(const HLTGenValObject& obj, edm::Handle<trigger::TriggerEvent>& triggerEvent) {
0142
0143 auto keyRange = getTrigObjIndicesOfCollection(*triggerEvent, collectionName_, hltProcessName_);
0144 const trigger::TriggerObject* bestMatch = nullptr;
0145 float bestDR2 = dR2limit_;
0146 for (trigger::size_type key = keyRange.first; key < keyRange.second; ++key) {
0147 if (passFilterSelection(key, *triggerEvent)) {
0148 const trigger::TriggerObject objTrig = triggerEvent->getObjects().at(key);
0149 if (isEventLevelVariable_) {
0150 bestDR2 = 0;
0151 bestMatch = &objTrig;
0152 break;
0153 } else {
0154 float dR2 = reco::deltaR2(obj, objTrig);
0155 if (dR2 < bestDR2) {
0156 bestDR2 = dR2;
0157 bestMatch = &objTrig;
0158 }
0159 }
0160 }
0161 }
0162 HLTGenValObject objWithTrig(obj);
0163 if (bestMatch) {
0164 objWithTrig.setTrigObject(*bestMatch);
0165 for (auto& hist : hists_)
0166 hist->fill(objWithTrig);
0167 }
0168 }
0169
0170
0171 void HLTGenResHistColl::book1D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig) {
0172
0173 auto vsVar = histConfig.getParameter<std::string>("vsVar");
0174 auto vsVarFunc = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVar);
0175 auto resVar = histConfig.getParameter<std::string>("resVar");
0176 auto resVarFunc = hltdqm::getUnaryFuncFloat<HLTGenValObject>(resVar);
0177
0178
0179 auto resBinLowEdgesDouble = histConfig.getParameter<std::vector<double>>("resBinLowEdges");
0180 auto vsBinLowEdgesDouble = histConfig.getParameter<std::vector<double>>("vsBinLowEdges");
0181
0182
0183 std::vector<edm::ParameterSet> allCutsVector = histConfig.getParameter<std::vector<edm::ParameterSet>>("rangeCuts");
0184 VarRangeCutColl<HLTGenValObject> rangeCuts(allCutsVector);
0185
0186
0187 const std::string tag;
0188
0189
0190 if (!vsVarFunc) {
0191 throw cms::Exception("ConfigError") << " vsVar " << vsVar << " is giving null ptr (likely empty) in " << __FILE__
0192 << "," << __LINE__ << std::endl;
0193 }
0194
0195
0196 std::vector<float> vsBinLowEdges = convertVec<double, float>(vsBinLowEdgesDouble);
0197 std::vector<float> resBinLowEdges = convertVec<double, float>(resBinLowEdgesDouble);
0198
0199 std::string histNameInc = getHistName(resVar);
0200 std::string histTitleInc = collectionName_ + "to " + objType_ + " " + resVar;
0201 if (!tag.empty()) {
0202 histNameInc += separator_ + tag;
0203 histTitleInc += " " + tag;
0204 }
0205 std::string histNameVsBase = getHistName(resVar, vsVar);
0206 std::string histName2D = histNameVsBase + separator_ + "2D";
0207 std::string histNameProfile = histNameVsBase + separator_ + "profile";
0208
0209 std::string histTitle2D = collectionName_ + "to " + objType_ + " " + resVar + " vs " + vsVar;
0210 std::string histTitleProfile = histTitle2D;
0211 if (!tag.empty()) {
0212 histName2D += separator_ + tag;
0213 histTitle2D += " " + tag;
0214 histNameProfile += separator_ + tag;
0215 histTitleProfile += " " + tag;
0216 }
0217
0218 auto resInc =
0219 iBooker.book1D(histNameInc.c_str(), histTitleInc.c_str(), resBinLowEdges.size() - 1, resBinLowEdges.data());
0220
0221 auto res2D = iBooker.book2D(histName2D.c_str(),
0222 histTitle2D.c_str(),
0223 vsBinLowEdges.size() - 1,
0224 vsBinLowEdges.data(),
0225 resBinLowEdges.size() - 1,
0226 resBinLowEdges.data());
0227
0228
0229 auto resProf = iBooker.bookProfile(histNameProfile.c_str(),
0230 histTitleProfile.c_str(),
0231 vsBinLowEdgesDouble.size() - 1,
0232 vsBinLowEdgesDouble.data(),
0233 0.2,
0234 5);
0235
0236 std::unique_ptr<HLTGenValHist> hist;
0237
0238 hist = std::make_unique<HLTGenValHist2D>(
0239 res2D->getTH2F(), resProf->getTProfile(), vsVar, resVar, vsVarFunc, resVarFunc, rangeCuts);
0240 hists_.emplace_back(std::move(hist));
0241 hist = std::make_unique<HLTGenValHist1D>(resInc->getTH1(), resVar, resVarFunc, rangeCuts);
0242 hists_.emplace_back(std::move(hist));
0243 }
0244
0245
0246 void HLTGenResHistColl::book2D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig2D) {
0247
0248 auto vsVarX = histConfig2D.getParameter<std::string>("vsVarX");
0249 auto vsVarY = histConfig2D.getParameter<std::string>("vsVarY");
0250 auto vsVarFuncX = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarX);
0251 auto vsVarFuncY = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarY);
0252 auto binLowEdgesDoubleX = histConfig2D.getParameter<std::vector<double>>("binLowEdgesX");
0253 auto binLowEdgesDoubleY = histConfig2D.getParameter<std::vector<double>>("binLowEdgesY");
0254
0255
0256 if (!vsVarFuncX) {
0257 throw cms::Exception("ConfigError") << " vsVar " << vsVarX << " is giving null ptr (likely empty) in " << __FILE__
0258 << "," << __LINE__ << std::endl;
0259 }
0260 if (!vsVarFuncY) {
0261 throw cms::Exception("ConfigError") << " vsVar " << vsVarY << " is giving null ptr (likely empty) in " << __FILE__
0262 << "," << __LINE__ << std::endl;
0263 }
0264
0265
0266 std::vector<float> binLowEdgesX;
0267 std::vector<float> binLowEdgesY;
0268 binLowEdgesX.reserve(binLowEdgesDoubleX.size());
0269 binLowEdgesY.reserve(binLowEdgesDoubleY.size());
0270 for (double lowEdge : binLowEdgesDoubleX)
0271 binLowEdgesX.push_back(lowEdge);
0272 for (double lowEdge : binLowEdgesDoubleY)
0273 binLowEdgesY.push_back(lowEdge);
0274
0275
0276 std::string histName =
0277 objType_ + separator_ + separator_ + "GEN" + separator_ + "2Dvs" + vsVarX + separator_ + vsVarY;
0278 std::string histTitle = objType_ + " GEN 2D vs " + vsVarX + " " + vsVarY;
0279
0280 auto me = iBooker.book2D(histName.c_str(),
0281 histTitle.c_str(),
0282 binLowEdgesX.size() - 1,
0283 binLowEdgesX.data(),
0284 binLowEdgesY.size() - 1,
0285 binLowEdgesY.data());
0286
0287 std::unique_ptr<HLTGenValHist> hist;
0288
0289 hist = std::make_unique<HLTGenValHist2D>(me->getTH2F(), vsVarX, vsVarY, vsVarFuncX, vsVarFuncY);
0290
0291 hists_.emplace_back(std::move(hist));
0292 }
0293
0294 std::string HLTGenResHistColl::getHistName(const std::string& resVar, const std::string& vsVar) const {
0295 std::string histName = histNamePrefix_ + separator_ + collectionName_ + separator_ + objType_ + separator_ + resVar;
0296 if (!vsVar.empty()) {
0297 histName += separator_ + "vs" + separator_ + vsVar;
0298 }
0299 return histName;
0300 }