Belle II Software  release-06-01-15
DQMHistAnalysisHLTModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <dqm/analysis/modules/DQMHistAnalysisHLTModule.h>
10 #include <framework/core/ModuleParam.templateDetails.h>
11 #include <TROOT.h>
12 
13 #include <hlt/utilities/Units.h>
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 namespace {
19  bool hasValue(const std::string& name, TH1* histogram)
20  {
21  return histogram->GetXaxis()->FindFixBin(name.c_str()) != -1;
22  }
23 
24  double getValue(const std::string& name, TH1* histogram)
25  {
26  if (not hasValue(name, histogram)) {
27  B2ERROR("This histogram does not have this value! (fallback value = -1)"
28  << LogVar("histogram", histogram->GetName())
29  << LogVar("value", name));
30  return -1;
31  }
32  auto binNumber = histogram->GetXaxis()->FindFixBin(name.c_str());
33  return histogram->GetBinContent(binNumber);
34  }
35 }
36 
37 REG_MODULE(DQMHistAnalysisHLT)
38 
40 {
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",
46  m_retentionPerUnit);
47 }
48 
49 void DQMHistAnalysisHLTModule::initialize()
50 {
51  // this seems to be important, or strange things happen
52  gROOT->cd();
53 
54  m_hEfficiency = {
55  new TCanvas("HLT/Ratio"),
56  new TH1F("Ratio", "Retention of selected HLT skims", 1, 0, 0)
57  };
58  m_hEfficiencyTotal = {
59  new TCanvas("HLT/RatioTotal"),
60  new TH1F("RatioTotal", "Ratio of Tags to all events", 1, 0, 0)
61  };
62  m_hCrossSection = {
63  new TCanvas("HLT/CrossSection"),
64  new TH1F("CrossSection", "Cross Section of triggered Events", 1, 0, 0)
65  };
66  m_hRatios = {
67  new TCanvas("HLT/RatioToBahbha"),
68  new TH1F("RatioToBahbha", "Ratio to bhabha events", 1, 0, 0)
69  };
70 
71  for (const std::string& l1Name : m_l1Histograms) {
72  m_hl1Ratios.emplace(l1Name, std::make_pair(
73  new TCanvas(("HLT/" + l1Name + "RatioToL1").c_str()),
74  // + 1 for total result
75  new TH1F((l1Name + "RatioToL1").c_str(), ("HLT Fractions for L1 " + l1Name).c_str(), 1, 0, 0)
76  ));
77  }
78 
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)
85  ));
86  }
87 
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");
96  }
97 
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");
105  }
106 
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);
113  histogram->Draw();
114  }
115 
116  m_hMeanTime = {
117  new TCanvas("HLT/MeanTime"),
118  new TH1F("MeanTime", "Mean processing time", 1, 0, 0)
119  };
120 
121  m_hErrorFlagFraction = {
122  new TCanvas("HLT/ErrorFlagFraction"),
123  new TH1D("ErrorFlagFraction", "Fraction of events with Error Flags", 1, 0, 0)
124  };
125 
126  m_hFilteredFractionPerUnit = {
127  new TCanvas("HLT/FilteredFractionPerUnit"),
128  new TH1D("FilteredFractionPerUnit", "Fraction of events filtered per unit", 1, 0, 0)
129  };
130 
131  m_hMeanBudgetTimePerUnit = {
132  new TCanvas("HLT/MeanBudgetTimePerUnit"),
133  new TH1F("MeanBudgetTimePerUnit", "Mean budget time per unit", 1, 0, 0)
134  };
135 
136  m_hMeanProcessingTimePerUnit = {
137  new TCanvas("HLT/MeanProcessingTimePerUnit"),
138  new TH1F("MeanProcessingTimePerUnit", "Mean processing time per unit", 1, 0, 0)
139  };
140 
141  m_hMeanMemory = {
142  new TCanvas("HLT/MeanMemoryChange"),
143  new TH1F("MeanMemoryChange", "Mean memory change [MB]", 1, 0, 0)
144  };
145 
146 #ifdef _BELLE2_EPICS
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");
151  }
152 #endif
153 }
154 
155 
156 void DQMHistAnalysisHLTModule::beginRun()
157 {
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;
160  canvas->Clear();
161  }
162 
163  for (auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
164  auto* canvas = canvasAndHisto.first;
165  canvas->Clear();
166  }
167 
168  for (auto& nameAndcanvasAndHisto : m_hl1Ratios) {
169  auto* canvas = nameAndcanvasAndHisto.second.first;
170  canvas->Clear();
171  }
172 
173  for (auto& nameAndcanvasAndHisto : m_hRetentionPerUnit) {
174  auto* canvas = nameAndcanvasAndHisto.second.first;
175  canvas->Clear();
176  }
177 }
178 
179 void DQMHistAnalysisHLTModule::event()
180 {
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");
192 
193  if (not filterHistogram) {
194  B2ERROR("Can not find the filter histogram!");
195  return;
196  }
197  if (not skimHistogram) {
198  B2ERROR("Can not find the skim histogram!");
199  return;
200  }
201  if (not totalResultHistogram) {
202  B2ERROR("Can not find the total result histogram!");
203  return;
204  }
205  if (not hltUnitNumberHistogram) {
206  B2ERROR("Can not find the HLT unit number histogram!");
207  return;
208  }
209  if (not processesPerUnitHistogram) {
210  B2ERROR("Can not find the processes per unit histogram!");
211  return;
212  }
213  if (not meanTimeHistogram) {
214  B2ERROR("Can not find the mean processing time histogram!");
215  return;
216  }
217  if (not errorFlagHistogram) {
218  B2ERROR("Can not find the error flag histogram!");
219  return;
220  }
221  if (not hltUnitNumberHistogram_filtered) {
222  B2ERROR("Can not find the HLT unit number after filter histogram!");
223  return;
224  }
225  if (not fullTimeMeanPerUnitHistogram) {
226  B2ERROR("Can not find the HLT budget time per unit histogram!");
227  return;
228  }
229  if (not processingTimeMeanPerUnitHistogram) {
230  B2ERROR("Can not find the HLT processing time per unit histogram!");
231  return;
232  }
233  if (not meanMemoryHistogram) {
234  B2ERROR("Can not find the mean memory change histogram!");
235  return;
236  }
237 
238  m_hEfficiencyTotal.second->Reset();
239  m_hEfficiency.second->Reset();
240  m_hCrossSection.second->Reset();
241  m_hRatios.second->Reset();
242 
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();
248 
249 #ifdef _BELLE2_EPICS
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");
253  } else {
254  instLuminosity = 0;
255  }
256 #endif
257 
258  m_hEfficiencyTotal.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfAllEvents);
259  if (instLuminosity != 0) {
260  m_hCrossSection.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfAllEvents * instLuminosity);
261  }
262  m_hRatios.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfBhabhaEvents);
263 
264  for (const auto& columnMapping : m_columnMapping) {
265  const auto& from = columnMapping.first;
266  const auto& to = columnMapping.second;
267 
268  double value = 0;
269  if (hasValue(from, filterHistogram)) {
270  value = getValue(from, filterHistogram);
271  } else if (hasValue(from, skimHistogram)) {
272  value = getValue(from, skimHistogram);
273  } else {
274  B2ERROR("Can not find value " << from << ". Will not use it!");
275  continue;
276  }
277 
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);
282  }
283  m_hRatios.second->Fill(to.c_str(), value / numberOfBhabhaEvents);
284  }
285 
286  for (const std::string& l1Name : m_l1Histograms) {
287  auto* histogram = m_hl1Ratios.at(l1Name).second;
288  histogram->Reset();
289 
290  auto* l1Histogram = findHist("softwaretrigger/" + l1Name);
291  auto* l1TotalResultHistogram = findHist("softwaretrigger/l1_total_result");
292 
293  if (not l1Histogram or not l1TotalResultHistogram) {
294  B2ERROR("Can not find L1 histograms from softwaretrigger!");
295  continue;
296  }
297 
298  for (const auto& columnMapping : m_columnMapping) {
299  const auto& from = columnMapping.first;
300  const auto& to = columnMapping.second;
301 
302  if (not hasValue(from, l1Histogram)) {
303  B2ERROR("Can not find label " << from << " in l1 histogram " << l1Name);
304  continue;
305  }
306 
307  if (not hasValue(l1Name, l1TotalResultHistogram)) {
308  B2ERROR("Can not find label " << l1Name << " in l1 total result histogram");
309  continue;
310  }
311 
312  const double hltValueInL1Bin = getValue(from, l1Histogram);
313  const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
314 
315  histogram->Fill(to.c_str(), hltValueInL1Bin / l1TotalResult);
316  }
317 
318  // and add a total result bin
319  const auto from = "hlt_result";
320  const auto to = "hlt_result";
321 
322  if (not hasValue(from, l1Histogram)) {
323  B2ERROR("Can not find label " << from << " in l1 histogram " << l1Name);
324  continue;
325  }
326 
327  if (not hasValue(l1Name, l1TotalResultHistogram)) {
328  B2ERROR("Can not find label " << l1Name << " in l1 total result histogram");
329  continue;
330  }
331 
332  const double hltValueInL1Bin = getValue(from, l1Histogram);
333  const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
334 
335  histogram->Fill(to, hltValueInL1Bin / l1TotalResult);
336  }
337 
338  for (const std::string& filterLine : m_retentionPerUnit) {
339  auto* histogram = m_hRetentionPerUnit.at(filterLine).second;
340  histogram->Reset();
341 
342  auto* filterLineHistogram = findHist("softwaretrigger/" + filterLine + "_per_unit");
343 
344  if (not filterLineHistogram) {
345  B2ERROR("Can not find " << filterLineHistogram << "_per_event histograms from softwaretrigger!");
346  continue;
347  }
348 
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);
353  } else {
354  double filterLineUnitValue = filterLineHistogram->GetBinContent(i);
355  histogram->SetBinContent(i, filterLineUnitValue / totalUnitValue);
356  }
357  }
358  }
359 
360  m_hMeanTime.second = (TH1F*) meanTimeHistogram->Clone("MeanTime");
361  m_hMeanTime.second->Scale(1 / numberOfProcesses);
362 
363  m_hMeanMemory.second = (TH1F*) meanMemoryHistogram->Clone("MeanMemoryChange");
364  m_hMeanMemory.second->Scale(1 / numberOfProcesses);
365 
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");
369 
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");
373 
374  m_hMeanBudgetTimePerUnit.second = (TH1F*) fullTimeMeanPerUnitHistogram->Clone("MeanBudgetTimePerUnit");
375  m_hMeanBudgetTimePerUnit.second->Divide(fullTimeMeanPerUnitHistogram, processesPerUnitHistogram);
376 
377  m_hMeanProcessingTimePerUnit.second = (TH1F*) processingTimeMeanPerUnitHistogram->Clone("MeanProcessingTimePerUnit");
378  m_hMeanProcessingTimePerUnit.second->Divide(processingTimeMeanPerUnitHistogram, processesPerUnitHistogram);
379 
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;
383 
384  canvas->cd();
385  histogram->LabelsDeflate("X");
386  histogram->Draw("hist");
387  canvas->Modified();
388  canvas->Update();
389  }
390 
391  for (auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
392  auto* canvas = canvasAndHisto.first;
393  auto* histogram = canvasAndHisto.second;
394 
395  canvas->cd();
396  histogram->LabelsDeflate("X");
397  histogram->Draw("hist");
398  histogram->SetStats(false);
399  canvas->Modified();
400  canvas->Update();
401  }
402 
403  for (auto& nameAndCanvasAndHisto : m_hl1Ratios) {
404  auto* canvas = nameAndCanvasAndHisto.second.first;
405  auto* histogram = nameAndCanvasAndHisto.second.second;
406 
407  canvas->cd();
408  histogram->LabelsDeflate("X");
409  histogram->Draw("hist");
410  canvas->Modified();
411  canvas->Update();
412  }
413 
414  for (auto& nameAndCanvasAndHisto : m_hRetentionPerUnit) {
415  auto* canvas = nameAndCanvasAndHisto.second.first;
416  auto* histogram = nameAndCanvasAndHisto.second.second;
417 
418  canvas->cd();
419  histogram->Draw("hist");
420  canvas->Modified();
421  canvas->Update();
422  }
423 }
424 
425 void DQMHistAnalysisHLTModule::terminate()
426 {
427 #ifdef _BELLE2_EPICS
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");
431  }
432 #endif
433 }
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.
Definition: Module.h:650
Abstract base class for different kinds of events.