Belle II Software  release-06-00-14
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 #include <dqm/analysis/modules/DQMHistAnalysisHLTModule.h>
9 #include <framework/core/ModuleParam.templateDetails.h>
10 #include <TROOT.h>
11 
12 #include <hlt/softwaretrigger/modules/dqm/SoftwareTriggerHLTDQMModule.h>
13 
14 using namespace std;
15 using namespace Belle2;
16 
17 namespace {
18  bool hasValue(const std::string& name, TH1* histogram)
19  {
20  return histogram->GetXaxis()->FindFixBin(name.c_str()) != -1;
21  }
22 
23  double getValue(const std::string& name, TH1* histogram)
24  {
25  if (not hasValue(name, histogram)) {
26  B2ERROR("This histogram does not have this value! (fallback value = -1)"
27  << LogVar("histogram", histogram->GetName())
28  << LogVar("value", name));
29  return -1;
30  }
31  auto binNumber = histogram->GetXaxis()->FindFixBin(name.c_str());
32  return histogram->GetBinContent(binNumber);
33  }
34 }
35 
36 REG_MODULE(DQMHistAnalysisHLT)
37 
38 DQMHistAnalysisHLTModule::DQMHistAnalysisHLTModule()
39 {
40  addParam("pvPrefix", m_pvPrefix, "EPICS PV Name for the inst. luminosity", m_pvPrefix);
41  addParam("bhabhaName", m_bhabhaName, "Name of the bhabha trigger to do a ratio against", m_bhabhaName);
42  addParam("columnMapping", m_columnMapping, "Which columns to use for calculating ratios and cross sections", m_columnMapping);
43  addParam("l1Histograms", m_l1Histograms, "Which l1 histograms to show", m_l1Histograms);
44  addParam("retentionPerUnit", m_retentionPerUnit, "Which HLT filter lines to use for calculation retention rate per unit",
45  m_retentionPerUnit);
46 }
47 
48 void DQMHistAnalysisHLTModule::initialize()
49 {
50  // this seems to be important, or strange things happen
51  gROOT->cd();
52 
53  m_hEfficiency = {
54  new TCanvas("HLT/Ratio"),
55  new TH1F("Ratio", "Retention of selected HLT skims", 1, 0, 0)
56  };
57  m_hEfficiencyTotal = {
58  new TCanvas("HLT/RatioTotal"),
59  new TH1F("RatioTotal", "Ratio of Tags to all events", 1, 0, 0)
60  };
61  m_hCrossSection = {
62  new TCanvas("HLT/CrossSection"),
63  new TH1F("CrossSection", "Cross Section of triggered Events", 1, 0, 0)
64  };
65  m_hRatios = {
66  new TCanvas("HLT/RatioToBahbha"),
67  new TH1F("RatioToBahbha", "Ratio to bhabha events", 1, 0, 0)
68  };
69 
70  for (const std::string& l1Name : m_l1Histograms) {
71  m_hl1Ratios.emplace(l1Name, std::make_pair(
72  new TCanvas(("HLT/" + l1Name + "RatioToL1").c_str()),
73  // + 1 for total result
74  new TH1F((l1Name + "RatioToL1").c_str(), ("HLT Fractions for L1 " + l1Name).c_str(), 1, 0, 0)
75  ));
76  }
77 
78  for (const std::string& filterLine : m_retentionPerUnit) {
79  m_hRetentionPerUnit.emplace(filterLine, std::make_pair(
80  new TCanvas(("HLT/" + filterLine + "_RetentionPerUnit").c_str()),
81  new TH1F((filterLine + "_RetentionPerUnit").c_str(), ("Retention rate per unit of: " + filterLine).c_str(),
82  SoftwareTrigger::HLTUnit::max_hlt_units + 1, 0,
83  SoftwareTrigger::HLTUnit::max_hlt_units + 1)
84  ));
85  }
86 
87  for (auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios}) {
88  auto* histogram = canvasAndHisto.second;
89  histogram->SetDirectory(0);
90  histogram->SetOption("bar");
91  histogram->SetFillStyle(0);
92  histogram->SetMinimum(0);
93  histogram->SetStats(false);
94  histogram->Draw("hist");
95  }
96 
97  for (auto& nameAndcanvasAndHisto : m_hl1Ratios) {
98  auto* histogram = nameAndcanvasAndHisto.second.second;
99  histogram->SetDirectory(0);
100  histogram->SetOption("bar");
101  histogram->SetFillStyle(0);
102  histogram->SetStats(false);
103  histogram->Draw("hist");
104  }
105 
106  for (auto& nameAndcanvasAndHisto : m_hRetentionPerUnit) {
107  auto* histogram = nameAndcanvasAndHisto.second.second;
108  histogram->SetDirectory(0);
109  histogram->SetOption("histe");
110  histogram->SetMinimum(0);
111  histogram->SetStats(false);
112  histogram->Draw();
113  }
114 
115  m_hMeanTime = {
116  new TCanvas("HLT/MeanTime"),
117  new TH1F("MeanTime", "Mean processing time", 1, 0, 0)
118  };
119 
120  m_hErrorFlagFraction = {
121  new TCanvas("HLT/ErrorFlagFraction"),
122  new TH1D("ErrorFlagFraction", "Fraction of events with Error Flags", 1, 0, 0)
123  };
124 
125  m_hFilteredFractionPerUnit = {
126  new TCanvas("HLT/FilteredFractionPerUnit"),
127  new TH1D("FilteredFractionPerUnit", "Fraction of events filtered per unit", 1, 0, 0)
128  };
129 
130  m_hMeanBudgetTimePerUnit = {
131  new TCanvas("HLT/MeanBudgetTimePerUnit"),
132  new TH1F("MeanBudgetTimePerUnit", "Mean budget time per unit", 1, 0, 0)
133  };
134 
135  m_hMeanProcessingTimePerUnit = {
136  new TCanvas("HLT/MeanProcessingTimePerUnit"),
137  new TH1F("MeanProcessingTimePerUnit", "Mean processing time per unit", 1, 0, 0)
138  };
139 
140  m_hMeanMemory = {
141  new TCanvas("HLT/MeanMemoryChange"),
142  new TH1F("MeanMemoryChange", "Mean memory change [MB]", 1, 0, 0)
143  };
144 
145 #ifdef _BELLE2_EPICS
146  if (not m_pvPrefix.empty()) {
147  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
148  SEVCHK(ca_create_channel(m_pvPrefix.data(), NULL, NULL, 10, &m_epicschid), "ca_create_channel failure");
149  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
150  }
151 #endif
152 }
153 
154 
155 void DQMHistAnalysisHLTModule::beginRun()
156 {
157  for (auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios, m_hMeanTime, m_hMeanBudgetTimePerUnit, m_hMeanProcessingTimePerUnit, m_hMeanMemory}) {
158  auto* canvas = canvasAndHisto.first;
159  canvas->Clear();
160  }
161 
162  for (auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
163  auto* canvas = canvasAndHisto.first;
164  canvas->Clear();
165  }
166 
167  for (auto& nameAndcanvasAndHisto : m_hl1Ratios) {
168  auto* canvas = nameAndcanvasAndHisto.second.first;
169  canvas->Clear();
170  }
171 
172  for (auto& nameAndcanvasAndHisto : m_hRetentionPerUnit) {
173  auto* canvas = nameAndcanvasAndHisto.second.first;
174  canvas->Clear();
175  }
176 }
177 
178 void DQMHistAnalysisHLTModule::event()
179 {
180  auto* filterHistogram = findHist("softwaretrigger/filter");
181  auto* skimHistogram = findHist("softwaretrigger/skim");
182  auto* totalResultHistogram = findHist("softwaretrigger/total_result");
183  auto* hltUnitNumberHistogram = findHist("softwaretrigger_before_filter/hlt_unit_number");
184  auto* processesPerUnitHistogram = findHist("timing_statistics/processesPerUnitHistogram");
185  auto* meanTimeHistogram = findHist("timing_statistics/meanTimeHistogram");
186  auto* errorFlagHistogram = findHist("softwaretrigger_before_filter/error_flag");
187  auto* hltUnitNumberHistogram_filtered = findHist("softwaretrigger/hlt_unit_number_after_filter");
188  auto* fullTimeMeanPerUnitHistogram = findHist("timing_statistics/fullTimeMeanPerUnitHistogram");
189  auto* processingTimeMeanPerUnitHistogram = findHist("timing_statistics/processingTimeMeanPerUnitHistogram");
190  auto* meanMemoryHistogram = findHist("timing_statistics/meanMemoryHistogram");
191 
192  if (not filterHistogram) {
193  B2ERROR("Can not find the filter histogram!");
194  return;
195  }
196  if (not skimHistogram) {
197  B2ERROR("Can not find the skim histogram!");
198  return;
199  }
200  if (not totalResultHistogram) {
201  B2ERROR("Can not find the total result histogram!");
202  return;
203  }
204  if (not hltUnitNumberHistogram) {
205  B2ERROR("Can not find the HLT unit number histogram!");
206  return;
207  }
208  if (not processesPerUnitHistogram) {
209  B2ERROR("Can not find the processes per unit histogram!");
210  return;
211  }
212  if (not meanTimeHistogram) {
213  B2ERROR("Can not find the mean processing time histogram!");
214  return;
215  }
216  if (not errorFlagHistogram) {
217  B2ERROR("Can not find the error flag histogram!");
218  return;
219  }
220  if (not hltUnitNumberHistogram_filtered) {
221  B2ERROR("Can not find the HLT unit number after filter histogram!");
222  return;
223  }
224  if (not fullTimeMeanPerUnitHistogram) {
225  B2ERROR("Can not find the HLT budget time per unit histogram!");
226  return;
227  }
228  if (not processingTimeMeanPerUnitHistogram) {
229  B2ERROR("Can not find the HLT processing time per unit histogram!");
230  return;
231  }
232  if (not meanMemoryHistogram) {
233  B2ERROR("Can not find the mean memory change histogram!");
234  return;
235  }
236 
237  m_hEfficiencyTotal.second->Reset();
238  m_hEfficiency.second->Reset();
239  m_hCrossSection.second->Reset();
240  m_hRatios.second->Reset();
241 
242  double instLuminosity = 0;
243  double numberOfAcceptedHLTEvents = getValue("total_result", totalResultHistogram);
244  double numberOfBhabhaEvents = getValue(m_bhabhaName, skimHistogram);
245  double numberOfAllEvents = hltUnitNumberHistogram->GetEntries();
246  double numberOfProcesses = processesPerUnitHistogram->GetEntries();
247 
248 #ifdef _BELLE2_EPICS
249  if (not m_pvPrefix.empty()) {
250  SEVCHK(ca_get(DBR_DOUBLE, m_epicschid, (void*)&instLuminosity), "ca_get failure");
251  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
252  } else {
253  instLuminosity = 0;
254  }
255 #endif
256 
257  m_hEfficiencyTotal.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfAllEvents);
258  if (instLuminosity != 0) {
259  m_hCrossSection.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfAllEvents * instLuminosity);
260  }
261  m_hRatios.second->Fill("total_result", numberOfAcceptedHLTEvents / numberOfBhabhaEvents);
262 
263  for (const auto& columnMapping : m_columnMapping) {
264  const auto& from = columnMapping.first;
265  const auto& to = columnMapping.second;
266 
267  double value = 0;
268  if (hasValue(from, filterHistogram)) {
269  value = getValue(from, filterHistogram);
270  } else if (hasValue(from, skimHistogram)) {
271  value = getValue(from, skimHistogram);
272  } else {
273  B2ERROR("Can not find value " << from << ". Will not use it!");
274  continue;
275  }
276 
277  m_hEfficiency.second->Fill(to.c_str(), value / numberOfAcceptedHLTEvents);
278  m_hEfficiencyTotal.second->Fill(to.c_str(), value / numberOfAllEvents);
279  if (instLuminosity != 0) {
280  m_hCrossSection.second->Fill(to.c_str(), value / numberOfAllEvents * instLuminosity);
281  }
282  m_hRatios.second->Fill(to.c_str(), value / numberOfBhabhaEvents);
283  }
284 
285  for (const std::string& l1Name : m_l1Histograms) {
286  auto* histogram = m_hl1Ratios.at(l1Name).second;
287  histogram->Reset();
288 
289  auto* l1Histogram = findHist("softwaretrigger/" + l1Name);
290  auto* l1TotalResultHistogram = findHist("softwaretrigger/l1_total_result");
291 
292  if (not l1Histogram or not l1TotalResultHistogram) {
293  B2ERROR("Can not find L1 histograms from softwaretrigger!");
294  continue;
295  }
296 
297  for (const auto& columnMapping : m_columnMapping) {
298  const auto& from = columnMapping.first;
299  const auto& to = columnMapping.second;
300 
301  if (not hasValue(from, l1Histogram)) {
302  B2ERROR("Can not find label " << from << " in l1 histogram " << l1Name);
303  continue;
304  }
305 
306  if (not hasValue(l1Name, l1TotalResultHistogram)) {
307  B2ERROR("Can not find label " << l1Name << " in l1 total result histogram");
308  continue;
309  }
310 
311  const double hltValueInL1Bin = getValue(from, l1Histogram);
312  const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
313 
314  histogram->Fill(to.c_str(), hltValueInL1Bin / l1TotalResult);
315  }
316 
317  // and add a total result bin
318  const auto from = "hlt_result";
319  const auto to = "hlt_result";
320 
321  if (not hasValue(from, l1Histogram)) {
322  B2ERROR("Can not find label " << from << " in l1 histogram " << l1Name);
323  continue;
324  }
325 
326  if (not hasValue(l1Name, l1TotalResultHistogram)) {
327  B2ERROR("Can not find label " << l1Name << " in l1 total result histogram");
328  continue;
329  }
330 
331  const double hltValueInL1Bin = getValue(from, l1Histogram);
332  const double l1TotalResult = getValue(l1Name, l1TotalResultHistogram);
333 
334  histogram->Fill(to, hltValueInL1Bin / l1TotalResult);
335  }
336 
337  for (const std::string& filterLine : m_retentionPerUnit) {
338  auto* histogram = m_hRetentionPerUnit.at(filterLine).second;
339  histogram->Reset();
340 
341  auto* filterLineHistogram = findHist("softwaretrigger/" + filterLine + "_per_unit");
342 
343  if (not filterLineHistogram) {
344  B2ERROR("Can not find " << filterLineHistogram << "_per_event histograms from softwaretrigger!");
345  continue;
346  }
347 
348  for (unsigned int i = 1; i <= SoftwareTrigger::HLTUnit::max_hlt_units + 1; i++) {
349  double totalUnitValue = hltUnitNumberHistogram->GetBinContent(i);
350  if (totalUnitValue == 0) {
351  histogram->Fill(i, 0);
352  } else {
353  double filterLineUnitValue = filterLineHistogram->GetBinContent(i);
354  histogram->SetBinContent(i, filterLineUnitValue / totalUnitValue);
355  }
356  }
357  }
358 
359  m_hMeanTime.second = (TH1F*) meanTimeHistogram->Clone("MeanTime");
360  m_hMeanTime.second->Scale(1 / numberOfProcesses);
361 
362  m_hMeanMemory.second = (TH1F*) meanMemoryHistogram->Clone("MeanMemoryChange");
363  m_hMeanMemory.second->Scale(1 / numberOfProcesses);
364 
365  m_hErrorFlagFraction.second = (TH1D*) errorFlagHistogram->Clone("ErrorFlagFraction");
366  m_hErrorFlagFraction.second->Scale(1 / numberOfAllEvents);
367  m_hErrorFlagFraction.second->SetTitle("Fraction of events with error flags");
368 
369  m_hFilteredFractionPerUnit.second = (TH1D*) hltUnitNumberHistogram_filtered->Clone("FilteredFractionPerUnit");
370  m_hFilteredFractionPerUnit.second->Divide(hltUnitNumberHistogram_filtered, hltUnitNumberHistogram);
371  m_hFilteredFractionPerUnit.second->SetTitle("Fraction of events filtered per unit");
372 
373  m_hMeanBudgetTimePerUnit.second = (TH1F*) fullTimeMeanPerUnitHistogram->Clone("MeanBudgetTimePerUnit");
374  m_hMeanBudgetTimePerUnit.second->Divide(fullTimeMeanPerUnitHistogram, processesPerUnitHistogram);
375 
376  m_hMeanProcessingTimePerUnit.second = (TH1F*) processingTimeMeanPerUnitHistogram->Clone("MeanProcessingTimePerUnit");
377  m_hMeanProcessingTimePerUnit.second->Divide(processingTimeMeanPerUnitHistogram, processesPerUnitHistogram);
378 
379  for (auto& canvasAndHisto : {m_hEfficiencyTotal, m_hEfficiency, m_hCrossSection, m_hRatios, m_hMeanTime, m_hMeanBudgetTimePerUnit, m_hMeanProcessingTimePerUnit, m_hMeanMemory}) {
380  auto* canvas = canvasAndHisto.first;
381  auto* histogram = canvasAndHisto.second;
382 
383  canvas->cd();
384  histogram->LabelsDeflate("X");
385  histogram->Draw("hist");
386  canvas->Modified();
387  canvas->Update();
388  }
389 
390  for (auto& canvasAndHisto : {m_hErrorFlagFraction, m_hFilteredFractionPerUnit}) {
391  auto* canvas = canvasAndHisto.first;
392  auto* histogram = canvasAndHisto.second;
393 
394  canvas->cd();
395  histogram->LabelsDeflate("X");
396  histogram->Draw("hist");
397  histogram->SetStats(false);
398  canvas->Modified();
399  canvas->Update();
400  }
401 
402  for (auto& nameAndCanvasAndHisto : m_hl1Ratios) {
403  auto* canvas = nameAndCanvasAndHisto.second.first;
404  auto* histogram = nameAndCanvasAndHisto.second.second;
405 
406  canvas->cd();
407  histogram->LabelsDeflate("X");
408  histogram->Draw("hist");
409  canvas->Modified();
410  canvas->Update();
411  }
412 
413  for (auto& nameAndCanvasAndHisto : m_hRetentionPerUnit) {
414  auto* canvas = nameAndCanvasAndHisto.second.first;
415  auto* histogram = nameAndCanvasAndHisto.second.second;
416 
417  canvas->cd();
418  histogram->Draw("hist");
419  canvas->Modified();
420  canvas->Update();
421  }
422 }
423 
424 void DQMHistAnalysisHLTModule::terminate()
425 {
426 #ifdef _BELLE2_EPICS
427  if (not m_pvPrefix.empty()) {
428  SEVCHK(ca_clear_channel(m_epicschid), "ca_clear_channel failure");
429  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
430  }
431 #endif
432 }
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.