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