Belle II Software  release-06-02-00
DQMHistDeltaHisto.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 // File : DQMHistDeltaHisto.cc
10 // Description : DQM Histogram analysis module to generate delta histograms
11 //-
12 
13 
14 #include <framework/core/ModuleParam.templateDetails.h>
15 #include <dqm/analysis/modules/DQMHistDeltaHisto.h>
16 #include <daq/slc/base/StringUtil.h>
17 #include <TROOT.h>
18 #include <TClass.h>
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 //-----------------------------------------------------------------
24 // Register the Module
25 //-----------------------------------------------------------------
26 REG_MODULE(DQMHistDeltaHisto)
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
34 {
35  addParam("Interval", m_interval, "Interval time for diff histos [s]", 180);
36  addParam("MonitoredHistos", m_monitored_histos, "List of histograms to monitor", vector<string>());
37  B2DEBUG(20, "DQMHistDeltaHisto: Constructor done.");
38 }
39 
40 
41 DQMHistDeltaHistoModule::~DQMHistDeltaHistoModule() { }
42 
43 void DQMHistDeltaHistoModule::initialize()
44 {
45  gROOT->cd();
46  B2DEBUG(20, "DQMHistDeltaHisto: initialized.");
47  for (auto& histoname : m_monitored_histos) {
48  queue<SSNODE*> hq;
49  m_histos_queues[histoname] = hq;
50  }
51  m_evtMetaDataPtr.isRequired();
52 }
53 
54 
55 void DQMHistDeltaHistoModule::beginRun()
56 {
57  B2DEBUG(20, "DQMHistDeltaHisto: beginRun called.");
58  for (auto& histoname : m_monitored_histos) {
59  queue<SSNODE*>& hq = m_histos_queues[histoname];
60  while (!hq.empty()) {
61  SSNODE* nn = hq.front();
62  clear_node(nn);
63  delete nn;
64  hq.pop();
65  }
66  }
67 }
68 
69 void DQMHistDeltaHistoModule::clear_node(SSNODE* n)
70 {
71  delete n->histo;
72  delete n->diff_histo;
73 }
74 
75 void DQMHistDeltaHistoModule::event()
76 {
77 
78  B2DEBUG(20, "DQMHistDeltaHisto: event called.");
79  if (!m_evtMetaDataPtr.isValid()) {
80  B2ERROR("No valid EventMetaData.");
81  return;
82  }
83  time_t cur_mtime = m_evtMetaDataPtr->getTime();
84 
85  for (auto& histoname : m_monitored_histos) {
86  TH1* hh = findHist(histoname.c_str());
87  if (hh == nullptr) continue;
88  if (hh->GetDimension() != 1) continue;
89  queue<SSNODE*>& hq = m_histos_queues[histoname];
90  if (hq.empty()) {
91  SSNODE* n = new SSNODE;
92  n->histo = (TH1*)hh->Clone();
93  n->diff_histo = (TH1*)hh->Clone();
94  n->time_modified = cur_mtime;
95  hq.push(n);
96  } else {
97  while (!hq.empty()) {
98  SSNODE* nn = hq.front();
99  if ((cur_mtime - nn->time_modified < m_interval) || (hq.size() == 1)) {
100  if (hq.back()->time_modified == cur_mtime) {
101  break;
102  }
103  SSNODE* n = new SSNODE;
104  n->histo = (TH1*)hh->Clone();
105  n->diff_histo = (TH1*)hh->Clone();
106  n->diff_histo->Add(nn->histo, -1);
107  n->time_modified = cur_mtime;
108  hq.push(n);
109  break;
110  } else {
111  clear_node(nn);
112  delete nn;
113  hq.pop();
114  }
115  }
116  }
117  TString a = histoname;
118  StringList s = StringUtil::split(a.Data(), '/');
119  std::string dirname = s[0];
120  std::string hname = s[1];
121  std::string canvas_name = dirname + "/c_" + hname;
122  TCanvas* c = findCanvas(canvas_name);
123  if (c == nullptr) continue;
124  TH1* h_diff = hq.back()->diff_histo;
125  h_diff->SetName((a + "_diff").Data());
126  if (h_diff->Integral() != 0) h_diff->Scale(hh->Integral() / h_diff->Integral());
127  c->cd();
128  h_diff->SetLineColor(kRed);
129  h_diff->SetLineStyle(kDotted);
130  h_diff->SetStats(kFALSE);
131  h_diff->Draw("hist,same");
132  c->Modified();
133  c->Update();
134  }
135 
136 }
137 
138 void DQMHistDeltaHistoModule::endRun()
139 {
140  B2DEBUG(20, "DQMHistDeltaHisto: endRun called");
141 }
142 
143 
144 void DQMHistDeltaHistoModule::terminate()
145 {
146  B2DEBUG(20, "DQMHistDeltaHisto: terminate called");
147 }
The base class for the histogram analysis module.
Class for generating snapshots for histograms.
#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.
The struct for the snapshots.
TH1 * histo
The histogram for snapshot.
time_t time_modified
Whether the histogram is not updated for a long time.