9 #include <dqm/analysis/modules/DQMHistAnalysisHLTModule.h>
10 #include <framework/core/ModuleParam.templateDetails.h>
13 #include <hlt/utilities/Units.h>
19 bool hasValue(
const std::string& name, TH1* histogram)
21 return histogram->GetXaxis()->FindFixBin(name.c_str()) != -1;
24 double getValue(
const std::string& name, TH1* histogram)
26 if (not hasValue(name, histogram)) {
27 B2ERROR(
"This histogram does not have this value! (fallback value = -1)"
28 <<
LogVar(
"histogram", histogram->GetName())
32 auto binNumber = histogram->GetXaxis()->FindFixBin(name.c_str());
33 return histogram->GetBinContent(binNumber);
41 addParam(
"pvPrefix", m_pvPrefix,
"EPICS PV Name for the inst. luminosity", m_pvPrefix);
42 addParam(
"bhabhaName", m_bhabhaName,
"Name of the bhabha trigger to do a ratio against", m_bhabhaName);
43 addParam(
"columnMapping", m_columnMapping,
"Which columns to use for calculating ratios and cross sections", m_columnMapping);
44 addParam(
"l1Histograms", m_l1Histograms,
"Which l1 histograms to show", m_l1Histograms);
45 addParam(
"retentionPerUnit", m_retentionPerUnit,
"Which HLT filter lines to use for calculation retention rate per unit",
49 void DQMHistAnalysisHLTModule::initialize()
55 new TCanvas(
"HLT/Ratio"),
56 new TH1F(
"Ratio",
"Retention of selected HLT skims", 1, 0, 0)
58 m_hEfficiencyTotal = {
59 new TCanvas(
"HLT/RatioTotal"),
60 new TH1F(
"RatioTotal",
"Ratio of Tags to all events", 1, 0, 0)
63 new TCanvas(
"HLT/CrossSection"),
64 new TH1F(
"CrossSection",
"Cross Section of triggered Events", 1, 0, 0)
67 new TCanvas(
"HLT/RatioToBahbha"),
68 new TH1F(
"RatioToBahbha",
"Ratio to bhabha events", 1, 0, 0)
71 for (
const std::string& l1Name : m_l1Histograms) {
72 m_hl1Ratios.emplace(l1Name, std::make_pair(
73 new TCanvas((
"HLT/" + l1Name +
"RatioToL1").c_str()),
75 new TH1F((l1Name +
"RatioToL1").c_str(), (
"HLT Fractions for L1 " + l1Name).c_str(), 1, 0, 0)
79 for (
const std::string& filterLine : m_retentionPerUnit) {
80 m_hRetentionPerUnit.emplace(filterLine, std::make_pair(
81 new TCanvas((
"HLT/" + filterLine +
"_RetentionPerUnit").c_str()),
82 new TH1F((filterLine +
"_RetentionPerUnit").c_str(), (
"Retention rate per unit of: " + filterLine).c_str(),
83 HLTUnits::max_hlt_units + 1, 0,
84 HLTUnits::max_hlt_units + 1)
88 for (
auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios}) {
89 auto* histogram = canvasAndHisto.second;
90 histogram->SetDirectory(0);
91 histogram->SetOption(
"bar");
92 histogram->SetFillStyle(0);
93 histogram->SetMinimum(0);
94 histogram->SetStats(
false);
95 histogram->Draw(
"hist");
98 for (
auto& nameAndcanvasAndHisto : m_hl1Ratios) {
99 auto* histogram = nameAndcanvasAndHisto.second.second;
100 histogram->SetDirectory(0);
101 histogram->SetOption(
"bar");
102 histogram->SetFillStyle(0);
103 histogram->SetStats(
false);
104 histogram->Draw(
"hist");
107 for (
auto& nameAndcanvasAndHisto : m_hRetentionPerUnit) {
108 auto* histogram = nameAndcanvasAndHisto.second.second;
109 histogram->SetDirectory(0);
110 histogram->SetOption(
"histe");
111 histogram->SetMinimum(0);
112 histogram->SetStats(
false);
117 new TCanvas(
"HLT/MeanTime"),
118 new TH1F(
"MeanTime",
"Mean processing time", 1, 0, 0)
121 m_hErrorFlagFraction = {
122 new TCanvas(
"HLT/ErrorFlagFraction"),
123 new TH1D(
"ErrorFlagFraction",
"Fraction of events with Error Flags", 1, 0, 0)
126 m_hFilteredFractionPerUnit = {
127 new TCanvas(
"HLT/FilteredFractionPerUnit"),
128 new TH1D(
"FilteredFractionPerUnit",
"Fraction of events filtered per unit", 1, 0, 0)
131 m_hMeanBudgetTimePerUnit = {
132 new TCanvas(
"HLT/MeanBudgetTimePerUnit"),
133 new TH1F(
"MeanBudgetTimePerUnit",
"Mean budget time per unit", 1, 0, 0)
136 m_hMeanProcessingTimePerUnit = {
137 new TCanvas(
"HLT/MeanProcessingTimePerUnit"),
138 new TH1F(
"MeanProcessingTimePerUnit",
"Mean processing time per unit", 1, 0, 0)
142 new TCanvas(
"HLT/MeanMemoryChange"),
143 new TH1F(
"MeanMemoryChange",
"Mean memory change [MB]", 1, 0, 0)
147 if (not m_pvPrefix.empty()) {
148 if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback),
"ca_context_create");
149 SEVCHK(ca_create_channel(m_pvPrefix.data(), NULL, NULL, 10, &m_epicschid),
"ca_create_channel failure");
150 SEVCHK(ca_pend_io(5.0),
"ca_pend_io failure");
156 void DQMHistAnalysisHLTModule::beginRun()
158 for (
auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios, m_hMeanTime, m_hMeanBudgetTimePerUnit, m_hMeanProcessingTimePerUnit, m_hMeanMemory}) {
159 auto* canvas = canvasAndHisto.first;
163 for (
auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
164 auto* canvas = canvasAndHisto.first;
168 for (
auto& nameAndcanvasAndHisto : m_hl1Ratios) {
169 auto* canvas = nameAndcanvasAndHisto.second.first;
173 for (
auto& nameAndcanvasAndHisto : m_hRetentionPerUnit) {
174 auto* canvas = nameAndcanvasAndHisto.second.first;
179 void DQMHistAnalysisHLTModule::event()
181 auto* filterHistogram = findHist(
"softwaretrigger/filter");
182 auto* skimHistogram = findHist(
"softwaretrigger/skim");
183 auto* totalResultHistogram = findHist(
"softwaretrigger/total_result");
184 auto* hltUnitNumberHistogram = findHist(
"softwaretrigger_before_filter/hlt_unit_number");
185 auto* processesPerUnitHistogram = findHist(
"timing_statistics/processesPerUnitHistogram");
186 auto* meanTimeHistogram = findHist(
"timing_statistics/meanTimeHistogram");
187 auto* errorFlagHistogram = findHist(
"softwaretrigger_before_filter/error_flag");
188 auto* hltUnitNumberHistogram_filtered = findHist(
"softwaretrigger/hlt_unit_number_after_filter");
189 auto* fullTimeMeanPerUnitHistogram = findHist(
"timing_statistics/fullTimeMeanPerUnitHistogram");
190 auto* processingTimeMeanPerUnitHistogram = findHist(
"timing_statistics/processingTimeMeanPerUnitHistogram");
191 auto* meanMemoryHistogram = findHist(
"timing_statistics/meanMemoryHistogram");
193 if (not filterHistogram) {
194 B2ERROR(
"Can not find the filter histogram!");
197 if (not skimHistogram) {
198 B2ERROR(
"Can not find the skim histogram!");
201 if (not totalResultHistogram) {
202 B2ERROR(
"Can not find the total result histogram!");
205 if (not hltUnitNumberHistogram) {
206 B2ERROR(
"Can not find the HLT unit number histogram!");
209 if (not processesPerUnitHistogram) {
210 B2ERROR(
"Can not find the processes per unit histogram!");
213 if (not meanTimeHistogram) {
214 B2ERROR(
"Can not find the mean processing time histogram!");
217 if (not errorFlagHistogram) {
218 B2ERROR(
"Can not find the error flag histogram!");
221 if (not hltUnitNumberHistogram_filtered) {
222 B2ERROR(
"Can not find the HLT unit number after filter histogram!");
225 if (not fullTimeMeanPerUnitHistogram) {
226 B2ERROR(
"Can not find the HLT budget time per unit histogram!");
229 if (not processingTimeMeanPerUnitHistogram) {
230 B2ERROR(
"Can not find the HLT processing time per unit histogram!");
233 if (not meanMemoryHistogram) {
234 B2ERROR(
"Can not find the mean memory change histogram!");
238 m_hEfficiencyTotal.second->Reset();
239 m_hEfficiency.second->Reset();
240 m_hCrossSection.second->Reset();
241 m_hRatios.second->Reset();
243 double instLuminosity = 0;
244 double numberOfAcceptedHLTEvents = getValue(
"total_result", totalResultHistogram);
245 double numberOfBhabhaEvents = getValue(m_bhabhaName, skimHistogram);
246 double numberOfAllEvents = hltUnitNumberHistogram->GetEntries();
247 double numberOfProcesses = processesPerUnitHistogram->GetEntries();
250 if (not m_pvPrefix.empty()) {
251 SEVCHK(ca_get(DBR_DOUBLE, m_epicschid, (
void*)&instLuminosity),
"ca_get failure");
252 SEVCHK(ca_pend_io(5.0),
"ca_pend_io failure");
258 m_hEfficiencyTotal.second->Fill(
"total_result", numberOfAcceptedHLTEvents / numberOfAllEvents);
259 if (instLuminosity != 0) {
260 m_hCrossSection.second->Fill(
"total_result", numberOfAcceptedHLTEvents / numberOfAllEvents * instLuminosity);
262 m_hRatios.second->Fill(
"total_result", numberOfAcceptedHLTEvents / numberOfBhabhaEvents);
264 for (
const auto& columnMapping : m_columnMapping) {
265 const auto& from = columnMapping.first;
266 const auto& to = columnMapping.second;
269 if (hasValue(from, filterHistogram)) {
270 value = getValue(from, filterHistogram);
271 }
else if (hasValue(from, skimHistogram)) {
272 value = getValue(from, skimHistogram);
274 B2ERROR(
"Can not find value " << from <<
". Will not use it!");
278 m_hEfficiency.second->Fill(to.c_str(), value / numberOfAcceptedHLTEvents);
279 m_hEfficiencyTotal.second->Fill(to.c_str(), value / numberOfAllEvents);
280 if (instLuminosity != 0) {
281 m_hCrossSection.second->Fill(to.c_str(), value / numberOfAllEvents * instLuminosity);
283 m_hRatios.second->Fill(to.c_str(), value / numberOfBhabhaEvents);
286 for (
const std::string& l1Name : m_l1Histograms) {
287 auto* histogram = m_hl1Ratios.at(l1Name).second;
290 auto* l1Histogram = findHist(
"softwaretrigger/" + l1Name);
291 auto* l1TotalResultHistogram = findHist(
"softwaretrigger/l1_total_result");
293 if (not l1Histogram or not l1TotalResultHistogram) {
294 B2ERROR(
"Can not find L1 histograms from softwaretrigger!");
298 for (
const auto& columnMapping : m_columnMapping) {
299 const auto& from = columnMapping.first;
300 const auto& to = columnMapping.second;
302 if (not hasValue(from, l1Histogram)) {
303 B2ERROR(
"Can not find label " << from <<
" in l1 histogram " << l1Name);
307 if (not hasValue(l1Name, l1TotalResultHistogram)) {
308 B2ERROR(
"Can not find label " << l1Name <<
" in l1 total result histogram");
312 const double hltValueInL1Bin = getValue(from, l1Histogram);
313 const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
315 histogram->Fill(to.c_str(), hltValueInL1Bin / l1TotalResult);
319 const auto from =
"hlt_result";
320 const auto to =
"hlt_result";
322 if (not hasValue(from, l1Histogram)) {
323 B2ERROR(
"Can not find label " << from <<
" in l1 histogram " << l1Name);
327 if (not hasValue(l1Name, l1TotalResultHistogram)) {
328 B2ERROR(
"Can not find label " << l1Name <<
" in l1 total result histogram");
332 const double hltValueInL1Bin = getValue(from, l1Histogram);
333 const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
335 histogram->Fill(to, hltValueInL1Bin / l1TotalResult);
338 for (
const std::string& filterLine : m_retentionPerUnit) {
339 auto* histogram = m_hRetentionPerUnit.at(filterLine).second;
342 auto* filterLineHistogram = findHist(
"softwaretrigger/" + filterLine +
"_per_unit");
344 if (not filterLineHistogram) {
345 B2ERROR(
"Can not find " << filterLineHistogram <<
"_per_event histograms from softwaretrigger!");
349 for (
unsigned int i = 1; i <= HLTUnits::max_hlt_units + 1; i++) {
350 double totalUnitValue = hltUnitNumberHistogram->GetBinContent(i);
351 if (totalUnitValue == 0) {
352 histogram->Fill(i, 0);
354 double filterLineUnitValue = filterLineHistogram->GetBinContent(i);
355 histogram->SetBinContent(i, filterLineUnitValue / totalUnitValue);
360 m_hMeanTime.second = (TH1F*) meanTimeHistogram->Clone(
"MeanTime");
361 m_hMeanTime.second->Scale(1 / numberOfProcesses);
363 m_hMeanMemory.second = (TH1F*) meanMemoryHistogram->Clone(
"MeanMemoryChange");
364 m_hMeanMemory.second->Scale(1 / numberOfProcesses);
366 m_hErrorFlagFraction.second = (TH1D*) errorFlagHistogram->Clone(
"ErrorFlagFraction");
367 m_hErrorFlagFraction.second->Scale(1 / numberOfAllEvents);
368 m_hErrorFlagFraction.second->SetTitle(
"Fraction of events with error flags");
370 m_hFilteredFractionPerUnit.second = (TH1D*) hltUnitNumberHistogram_filtered->Clone(
"FilteredFractionPerUnit");
371 m_hFilteredFractionPerUnit.second->Divide(hltUnitNumberHistogram_filtered, hltUnitNumberHistogram);
372 m_hFilteredFractionPerUnit.second->SetTitle(
"Fraction of events filtered per unit");
374 m_hMeanBudgetTimePerUnit.second = (TH1F*) fullTimeMeanPerUnitHistogram->Clone(
"MeanBudgetTimePerUnit");
375 m_hMeanBudgetTimePerUnit.second->Divide(fullTimeMeanPerUnitHistogram, processesPerUnitHistogram);
377 m_hMeanProcessingTimePerUnit.second = (TH1F*) processingTimeMeanPerUnitHistogram->Clone(
"MeanProcessingTimePerUnit");
378 m_hMeanProcessingTimePerUnit.second->Divide(processingTimeMeanPerUnitHistogram, processesPerUnitHistogram);
380 for (
auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios, m_hMeanTime, m_hMeanBudgetTimePerUnit, m_hMeanProcessingTimePerUnit, m_hMeanMemory}) {
381 auto* canvas = canvasAndHisto.first;
382 auto* histogram = canvasAndHisto.second;
385 histogram->LabelsDeflate(
"X");
386 histogram->Draw(
"hist");
391 for (
auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
392 auto* canvas = canvasAndHisto.first;
393 auto* histogram = canvasAndHisto.second;
396 histogram->LabelsDeflate(
"X");
397 histogram->Draw(
"hist");
398 histogram->SetStats(
false);
403 for (
auto& nameAndCanvasAndHisto : m_hl1Ratios) {
404 auto* canvas = nameAndCanvasAndHisto.second.first;
405 auto* histogram = nameAndCanvasAndHisto.second.second;
408 histogram->LabelsDeflate(
"X");
409 histogram->Draw(
"hist");
414 for (
auto& nameAndCanvasAndHisto : m_hRetentionPerUnit) {
415 auto* canvas = nameAndCanvasAndHisto.second.first;
416 auto* histogram = nameAndCanvasAndHisto.second.second;
419 histogram->Draw(
"hist");
425 void DQMHistAnalysisHLTModule::terminate()
428 if (not m_pvPrefix.empty()) {
429 SEVCHK(ca_clear_channel(m_epicschid),
"ca_clear_channel failure");
430 SEVCHK(ca_pend_io(5.0),
"ca_pend_io failure");
Class for HLT-related histogram analysis.
Class to store variables with their name which were sent to the logging service.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.