File indexing completed on 2024-04-06 12:32:41
0001 #include "Validation/HLTrigger/interface/HLTGenValHistCollFilter.h"
0002
0003
0004 HLTGenValHistCollFilter::HLTGenValHistCollFilter(edm::ParameterSet filterCollConfig) {
0005 objType_ = filterCollConfig.getParameter<std::string>("objType");
0006 filter_ = filterCollConfig.getParameter<std::string>("filterName");
0007 path_ = filterCollConfig.getParameter<std::string>("pathName");
0008 hltProcessName_ = filterCollConfig.getParameter<std::string>("hltProcessName");
0009 dR2limit_ = filterCollConfig.getParameter<double>("dR2limit");
0010 separator_ = "__";
0011 }
0012
0013 edm::ParameterSetDescription HLTGenValHistCollFilter::makePSetDescription() {
0014 edm::ParameterSetDescription desc;
0015 desc.add<std::string>("objType", "");
0016 desc.add<std::string>("hltProcessName", "HLT");
0017 desc.add<double>("dR2limit", 0.1);
0018 desc.add<std::string>("filterName", "");
0019 desc.add<std::string>("pathName", "");
0020 return desc;
0021 }
0022
0023
0024 void HLTGenValHistCollFilter::bookHists(DQMStore::IBooker& iBooker,
0025 const std::vector<edm::ParameterSet>& histConfigs,
0026 const std::vector<edm::ParameterSet>& histConfigs2D) {
0027 for (const auto& histConfig : histConfigs)
0028 book1D(iBooker, histConfig);
0029 for (const auto& histConfig : histConfigs2D)
0030 book2D(iBooker, histConfig);
0031 }
0032
0033
0034 void HLTGenValHistCollFilter::fillHists(const HLTGenValObject& obj, edm::Handle<trigger::TriggerEvent>& triggerEvent) {
0035
0036
0037 if (filter_ == "beforeAnyFilter") {
0038 for (auto& hist : hists_)
0039 hist->fill(obj);
0040 } else {
0041
0042
0043
0044 edm::InputTag filterTag(filter_, "", hltProcessName_);
0045 size_t filterIndex = triggerEvent->filterIndex(filterTag);
0046
0047
0048 trigger::TriggerObjectCollection allTriggerObjects = triggerEvent->getObjects();
0049 trigger::TriggerObjectCollection selectedObjects;
0050 if (filterIndex < triggerEvent->sizeFilters()) {
0051 const auto& keys = triggerEvent->filterKeys(filterIndex);
0052 for (unsigned short key : keys) {
0053 trigger::TriggerObject foundObject = allTriggerObjects[key];
0054 selectedObjects.push_back(foundObject);
0055 }
0056 }
0057
0058
0059 const static std::vector<std::string> eventLevelVariables = {"AK4HT", "AK8HT", "MET"};
0060 if (std::find(eventLevelVariables.begin(), eventLevelVariables.end(), objType_) != eventLevelVariables.end()) {
0061
0062 if (!selectedObjects.empty())
0063 for (auto& hist : hists_)
0064 hist->fill(obj);
0065 } else {
0066
0067 double mindR2 = 99999;
0068 for (const auto& filterobj : selectedObjects) {
0069 double dR = deltaR2(obj, filterobj);
0070 if (dR < mindR2)
0071 mindR2 = dR;
0072 }
0073 if (mindR2 < dR2limit_)
0074 for (auto& hist : hists_)
0075 hist->fill(obj);
0076 }
0077 }
0078 }
0079
0080
0081 void HLTGenValHistCollFilter::book1D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig) {
0082
0083 auto vsVar = histConfig.getParameter<std::string>("vsVar");
0084 auto vsVarFunc = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVar);
0085
0086
0087 HLTGenValPathSpecificSettingParser parser =
0088 HLTGenValPathSpecificSettingParser(histConfig.getParameter<std::string>("pathSpecificCuts"),
0089 histConfig.getParameter<std::vector<edm::ParameterSet>>("binnings"),
0090 vsVar);
0091
0092
0093 auto binLowEdgesDouble = histConfig.getParameter<std::vector<double>>("binLowEdges");
0094 if (parser.havePathSpecificBins())
0095 binLowEdgesDouble = *parser.getPathSpecificBins();
0096
0097
0098 std::vector<edm::ParameterSet> allCutsVector = histConfig.getParameter<std::vector<edm::ParameterSet>>("rangeCuts");
0099 std::vector<edm::ParameterSet> pathSpecificCutsVector = *parser.getPathSpecificCuts();
0100 allCutsVector.insert(allCutsVector.end(), pathSpecificCutsVector.begin(), pathSpecificCutsVector.end());
0101 VarRangeCutColl<HLTGenValObject> rangeCuts(allCutsVector);
0102
0103
0104 const std::string tag = *parser.getTag();
0105
0106
0107 if (!vsVarFunc) {
0108 throw cms::Exception("ConfigError") << " vsVar " << vsVar << " is giving null ptr (likely empty) in " << __FILE__
0109 << "," << __LINE__ << std::endl;
0110 }
0111
0112
0113 std::vector<float> binLowEdges;
0114 binLowEdges.reserve(binLowEdgesDouble.size());
0115 for (double lowEdge : binLowEdgesDouble)
0116 binLowEdges.push_back(lowEdge);
0117
0118
0119
0120
0121 std::string filterName = filter_;
0122 if (filterName.rfind('-', 0) == 0)
0123 filterName.erase(0, 1);
0124
0125 std::string histName, histTitle;
0126 if (filter_ == "beforeAnyFilter") {
0127 histName = objType_ + separator_ + path_ + separator_ + "GEN" + separator_ + "vs" + vsVar;
0128 histTitle = objType_ + " " + path_ + " GEN vs " + vsVar;
0129 if (!tag.empty()) {
0130 histName += separator_ + tag;
0131 histTitle += " " + tag;
0132 }
0133 } else {
0134 histName = objType_ + separator_ + path_ + separator_ + filterName + separator_ + "vs" + vsVar;
0135 histTitle = objType_ + " " + path_ + " " + filterName + " vs" + vsVar;
0136
0137
0138 if (!tag.empty()) {
0139 histName += separator_ + tag;
0140 histTitle += " " + tag;
0141 }
0142 }
0143
0144 auto me = iBooker.book1D(
0145 histName.c_str(), histTitle.c_str(), binLowEdges.size() - 1, binLowEdges.data());
0146
0147 std::unique_ptr<HLTGenValHist> hist;
0148
0149 hist = std::make_unique<HLTGenValHist1D>(me->getTH1(), vsVar, vsVarFunc, rangeCuts);
0150
0151 hists_.emplace_back(std::move(hist));
0152 }
0153
0154
0155 void HLTGenValHistCollFilter::book2D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig2D) {
0156
0157 auto vsVarX = histConfig2D.getParameter<std::string>("vsVarX");
0158 auto vsVarY = histConfig2D.getParameter<std::string>("vsVarY");
0159 auto vsVarFuncX = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarX);
0160 auto vsVarFuncY = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarY);
0161 auto binLowEdgesDoubleX = histConfig2D.getParameter<std::vector<double>>("binLowEdgesX");
0162 auto binLowEdgesDoubleY = histConfig2D.getParameter<std::vector<double>>("binLowEdgesY");
0163
0164
0165 if (!vsVarFuncX) {
0166 throw cms::Exception("ConfigError") << " vsVar " << vsVarX << " is giving null ptr (likely empty) in " << __FILE__
0167 << "," << __LINE__ << std::endl;
0168 }
0169 if (!vsVarFuncY) {
0170 throw cms::Exception("ConfigError") << " vsVar " << vsVarY << " is giving null ptr (likely empty) in " << __FILE__
0171 << "," << __LINE__ << std::endl;
0172 }
0173
0174
0175 std::vector<float> binLowEdgesX;
0176 std::vector<float> binLowEdgesY;
0177 binLowEdgesX.reserve(binLowEdgesDoubleX.size());
0178 binLowEdgesY.reserve(binLowEdgesDoubleY.size());
0179 for (double lowEdge : binLowEdgesDoubleX)
0180 binLowEdgesX.push_back(lowEdge);
0181 for (double lowEdge : binLowEdgesDoubleY)
0182 binLowEdgesY.push_back(lowEdge);
0183
0184
0185
0186
0187 std::string filterName = filter_;
0188 if (filterName.rfind('-', 0) == 0)
0189 filterName.erase(0, 1);
0190
0191 std::string histName, histTitle;
0192 if (filter_ == "beforeAnyFilter") {
0193 histName = objType_ + separator_ + path_ + separator_ + "GEN" + separator_ + "2Dvs" + vsVarX + separator_ + vsVarY;
0194 histTitle = objType_ + " " + path_ + " GEN 2D vs " + vsVarX + " " + vsVarY;
0195 } else {
0196 histName = objType_ + separator_ + path_ + separator_ + filterName + separator_ + "2Dvs" + vsVarX + vsVarY;
0197 histTitle = objType_ + " " + path_ + " " + filterName + " 2D vs" + vsVarX + " " + vsVarY;
0198 }
0199
0200 auto me = iBooker.book2D(histName.c_str(),
0201 histTitle.c_str(),
0202 binLowEdgesX.size() - 1,
0203 binLowEdgesX.data(),
0204 binLowEdgesY.size() - 1,
0205 binLowEdgesY.data());
0206
0207 std::unique_ptr<HLTGenValHist> hist;
0208
0209 hist = std::make_unique<HLTGenValHist2D>(me->getTH2F(), vsVarX, vsVarY, vsVarFuncX, vsVarFuncY);
0210
0211 hists_.emplace_back(std::move(hist));
0212 }