Belle II Software  release-05-01-25
DQMHistSnapshots.cc
1 //+
2 // File : DQMHistSnapshots.cc
3 // Description : DQM Histogram analysis module, generate snapshots of histograms
4 //
5 // Author : Boqun Wang, U. of Cincinnati
6 // Date : yesterday
7 //-
8 
9 
10 #include <framework/core/ModuleParam.templateDetails.h>
11 #include <dqm/analysis/modules/DQMHistSnapshots.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(DQMHistSnapshots)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
32 {
33  addParam("CheckInterval", m_check_interval, "Interval between two checks [s]", 180);
34  B2DEBUG(1, "DQMHistSnapshots: Constructor done.");
35 }
36 
37 
38 DQMHistSnapshotsModule::~DQMHistSnapshotsModule() { }
39 
40 void DQMHistSnapshotsModule::initialize()
41 {
42  gROOT->cd();
43  B2DEBUG(20, "DQMHistSnapshots: initialized.");
44 }
45 
46 
47 void DQMHistSnapshotsModule::beginRun()
48 {
49  m_ssnode.clear();
50  B2DEBUG(20, "DQMHistSnapshots: beginRun called.");
51 }
52 
53 DQMHistSnapshotsModule::SSNODE* DQMHistSnapshotsModule::find_snapshot(TString a)
54 {
55  for (auto& it : m_ssnode) {
56  if (it->histo->GetName() == a)
57  return it;
58  }
59  return NULL;
60 }
61 
62 TCanvas* DQMHistSnapshotsModule::find_canvas(TString canvas_name)
63 {
64  TIter nextkey(gROOT->GetListOfCanvases());
65  TObject* obj = NULL;
66 
67  while ((obj = (TObject*)nextkey())) {
68  if (obj->IsA()->InheritsFrom("TCanvas")) {
69  if (obj->GetName() == canvas_name)
70  return (TCanvas*)obj;
71  }
72  }
73  return NULL;
74 }
75 
76 void DQMHistSnapshotsModule::event()
77 {
78 
79  time_t cur_time = time(NULL);
80  int check = 0;
81  if ((m_last_check == 0) || (cur_time - m_last_check > m_check_interval)) {
82  check = 1;
83  m_last_check = cur_time;
84  }
85 
86  const HistList& hlist = getHistList();
87 
88  for (HistList::const_iterator it = hlist.begin(); it != hlist.end(); ++it) {
89  TString a = it->first;
90 
91  SSNODE* n = find_snapshot(a);
92  if (n == NULL) { // no existing snapshot, create new one
93  n = new SSNODE;
94  n->histo = (TH1*) it->second->Clone();
95 
96  StringList s = StringUtil::split(a.Data(), '/');
97  std::string dirname = s[0];
98  std::string hname = s[1];
99  std::string canvas_name = dirname + "/c_" + hname;
100  n->canvas = find_canvas(canvas_name);
101  n->stale = 0;
102 
103  m_ssnode.push_back(n);
104  } else {
105  TH1* h = it->second;
106  if (check == 1) {
107  if (h->GetEntries() > n->histo->GetEntries()) { // histogram has been updated
108  delete n->histo;
109  n->histo = (TH1*)h->Clone();
110  n->stale = 0;
111  } else { // notify that the histogram is stale
112  n->stale = 1;
113  }
114  }
115  if (n->stale == 1 && n->canvas != NULL) {
116  h->SetTitle((h->GetTitle() + string(" [STALLED]")).c_str());
117  }
118  }
119 
120  }
121 }
122 
123 void DQMHistSnapshotsModule::endRun()
124 {
125  B2DEBUG(20, "DQMHistSnapshots: endRun called");
126 }
127 
128 
129 void DQMHistSnapshotsModule::terminate()
130 {
131  B2DEBUG(20, "DQMHistSnapshots: terminate called");
132 }
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistSnapshotsModule
Class for generating snapshots for histograms.
Definition: DQMHistSnapshots.h:31
Belle2::DQMHistSnapshotsModule::SSNODE
The struct for the snapshots.
Definition: DQMHistSnapshots.h:36
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisModule::HistList
std::map< std::string, TH1 * > HistList
The type of list of histograms.
Definition: DQMHistAnalysis.h:66
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27