Belle II Software  release-05-01-25
DQMHistAnalysisHLTMonObj.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Luka Santelj *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // Own include
12 #include <dqm/analysis/modules/DQMHistAnalysisHLTMonObj.h>
13 // software trigger include
14 #include <hlt/softwaretrigger/modules/dqm/SoftwareTriggerHLTDQMModule.h>
15 //DQM
16 #include <dqm/analysis/modules/DQMHistAnalysis.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 
21 //-----------------------------------------------------------------
22 // Register module
23 //-----------------------------------------------------------------
24 
25 REG_MODULE(DQMHistAnalysisHLTMonObj);
26 
27 DQMHistAnalysisHLTMonObjModule::DQMHistAnalysisHLTMonObjModule()
29 {
30  setDescription("Produces MonitoringObject for the HLT from the available DQM histograms");
32 }
33 
35 {
36 }
37 
39 {
40  // make monitoring object related to this module
41  // if monitoring object already exists this will return pointer to it
43 
44 }
45 
47 {
48 }
49 
51 {
52  // can put the analysis code here or in endRun() function
53  // for the start tests we will store output only end of run so better to put code there
54 }
55 
57 {
58 
59  // get existing histograms produced by DQM modules
60  TH1* h_hlt = findHist("softwaretrigger/total_result");
61  TH1* h_skim = findHist("softwaretrigger/skim");
62  TH1* h_budget = findHist("timing_statistics/fullTimeHistogram");
63  TH1* h_processing = findHist("timing_statistics/processingTimeHistogram");
64  TH1* h_l1 = findHist("softwaretrigger_before_filter/hlt_unit_number");
65 
66  double n_hlt = 0.;
67  if (h_hlt) n_hlt = (double)h_hlt->GetBinContent((h_hlt->GetXaxis())->FindFixBin("total_result"));
68  m_monObj->setVariable("n_hlt", n_hlt);
69  double n_l1 = 0.;
70  if (h_l1) n_l1 = h_l1->GetEntries();
71  m_monObj->setVariable("n_l1", n_l1);
72 
73  if (h_skim) {
74  // loop bins, add variable to monObj named as "effCS_" + bin label w/o "accept"
75  for (int ibin = 1; ibin < h_skim->GetXaxis()->GetNbins() + 1; ibin++) {
76  double nentr = (double)h_skim->GetBinContent(ibin);
77  std::string bin_name(h_skim->GetXaxis()->GetBinLabel(ibin));
78  m_monObj->setVariable(bin_name.replace(0, 6, "effCS"), nentr);
79  }
80  }
81 
82  double bgt = 0.;
83  if (h_budget) bgt = h_budget->GetMean();
84  m_monObj->setVariable("budget_time", bgt);
85 
86  m_monObj->setVariable("n_l1_x_budget_time", n_l1 * bgt);
87 
88  double procTime = 0.;
89  if (h_processing) procTime = h_processing->GetMean();
90  m_monObj->setVariable("processing_time", procTime);
91 
92  TH1* h_budgetUnit = nullptr;
93 
94  for (unsigned int index = 1; index <= SoftwareTrigger::HLTUnit::max_hlt_units; index++) {
95  // add budget time per unit
96  h_budgetUnit = findHist(("timing_statistics/fullTimePerUnitHistogram_HLT" + std::to_string(index)).c_str());
97  double bgunit = 0.;
98  if (h_budgetUnit) bgunit = h_budgetUnit->GetMean();
99  m_monObj->setVariable(("budget_time_HLT" + std::to_string(index)).c_str(), bgunit);
100  // add processing time per unit
101  h_budgetUnit = findHist(("timing_statistics/processingTimePerUnitHistogram_HLT" + std::to_string(index)).c_str());
102  if (h_budgetUnit) bgunit = h_budgetUnit->GetMean();
103  else bgunit = 0.;
104  m_monObj->setVariable(("processing_time_HLT" + std::to_string(index)).c_str(), bgunit);
105  }
106 
107  B2DEBUG(20, "DQMHistAnalysisHLTMonObj : endRun called");
108 }
109 
111 {
112  B2DEBUG(20, "terminate called");
113 }
Belle2::SoftwareTrigger::HLTUnit::max_hlt_units
static constexpr unsigned int max_hlt_units
Maximum number of HLT units used during the experiment.
Definition: SoftwareTriggerHLTDQMModule.h:135
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::Module::c_ParallelProcessingCertified
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:82
Belle2::DQMHistAnalysisHLTMonObjModule::initialize
virtual void initialize() override
Initialize the Module.
Definition: DQMHistAnalysisHLTMonObj.cc:38
Belle2::DQMHistAnalysisModule::findHist
static TH1 * findHist(const std::string &histname)
Find histogram.
Definition: DQMHistAnalysis.cc:83
Belle2::DQMHistAnalysisHLTMonObjModule::~DQMHistAnalysisHLTMonObjModule
virtual ~DQMHistAnalysisHLTMonObjModule()
Destructor.
Definition: DQMHistAnalysisHLTMonObj.cc:34
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::MonitoringObject::setVariable
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
Definition: MonitoringObject.h:134
Belle2::DQMHistAnalysisHLTMonObjModule::m_monObj
MonitoringObject * m_monObj
MonitoringObject to be produced by this module.
Definition: DQMHistAnalysisHLTMonObj.h:74
Belle2::DQMHistAnalysisHLTMonObjModule::terminate
virtual void terminate() override
Termination action.
Definition: DQMHistAnalysisHLTMonObj.cc:110
Belle2::DQMHistAnalysisHLTMonObjModule::endRun
virtual void endRun() override
End-of-run action.
Definition: DQMHistAnalysisHLTMonObj.cc:56
Belle2::DQMHistAnalysisHLTMonObjModule::event
virtual void event() override
Event processor.
Definition: DQMHistAnalysisHLTMonObj.cc:50
Belle2::DQMHistAnalysisModule::getMonitoringObject
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
Definition: DQMHistAnalysis.cc:55
Belle2::DQMHistAnalysisHLTMonObjModule::beginRun
virtual void beginRun() override
Begin run function.
Definition: DQMHistAnalysisHLTMonObj.cc:46
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27