Belle II Software release-09-00-00
DQMHistReferenceModule.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#include <dqm/analysis/modules/DQMHistReferenceModule.h>
10#include <TROOT.h>
11#include <TStyle.h>
12#include <TClass.h>
13#include <TKey.h>
14
15using namespace std;
16using namespace Belle2;
17
18//-----------------------------------------------------------------
19// Register the Module
20//-----------------------------------------------------------------
21REG_MODULE(DQMHistReference);
22
23//-----------------------------------------------------------------
24// Implementation
25//-----------------------------------------------------------------
26
28{
29 //Parameter definition
30 addParam("ReferenceFile", m_referenceFile, "Name of the reference histrogram files", string(""));
31 B2DEBUG(1, "DQMHistReference: Constructor done.");
32}
33
34
36
38{
39 gStyle->SetOptStat(0);
40 gStyle->SetStatStyle(1);
41 gStyle->SetOptDate(22);// Date and Time in Bottom Right, does no work
42
43 B2DEBUG(1, "DQMHistReference: initialized.");
44}
45
47{
48 B2DEBUG(1, "DQMHistReference: beginRun called.");
49 m_firstInRun = true;
50}
51
53{
54 B2DEBUG(1, "DQMHistReference: reading references from input root file");
55
56 string run_type = getRunType();
57 if (run_type == "") run_type = "default";
58
59 B2INFO("DQMHistReference: run_type " << run_type);
60
61 for (auto& it : m_pnode) {
62 // clear ref histos from memory
63 if (it.ref_org) delete it.ref_org;
64 if (it.ref_clone) delete it.ref_clone;
65 }
66 m_pnode.clear();
67 B2INFO("DQMHistReference: clear m_pnode. size: " << m_pnode.size());
68
69 TFile* refFile = new TFile(m_referenceFile.c_str(), "READ");
70
71 if (refFile->IsZombie()) {
72 B2INFO("DQMHistReference: reference file " << m_referenceFile << " does not exist. No references will be used!");
73 refFile->Close();
74 delete refFile;
75 return;
76 }
77
78 B2INFO("DQMHistReference: use reference file " << m_referenceFile);
79
80 TIter nextkey(refFile->GetListOfKeys());
81 TKey* key;
82 while ((key = (TKey*)nextkey())) {
83 if (key->IsFolder() && string(key->GetName()) == string("ref")) {
84 TDirectory* refdir = (TDirectory*)key->ReadObj(); // ReadObj -> I own it
85 TIter nextDetDir(refdir->GetListOfKeys());
86 TKey* detDir;
87 // detector folders
88 while ((detDir = (TKey*)nextDetDir())) {
89 if (!detDir->IsFolder()) continue;
90 TIter nextTypeDir(((TDirectory*)detDir->ReadObj())->GetListOfKeys());
91 TKey* typeDir;
92 TDirectory* foundDir = NULL;
93 // run type folders (get the run type corresponding folder or use default one)
94 while ((typeDir = (TKey*)nextTypeDir())) {
95 if (!typeDir->IsFolder()) continue;
96 if (string(typeDir->GetName()) == run_type) {
97 foundDir = (TDirectory*)typeDir->ReadObj(); // ReadObj -> I own it
98 break;
99 }
100 if (string(typeDir->GetName()) == "default") foundDir = (TDirectory*)typeDir->ReadObj(); // ReadObj -> I own it
101 }
102 string dirname = detDir->GetName();
103 if (!foundDir) {
104 B2INFO("No run type specific or default references available for " << dirname);
105 } else {
106 B2INFO("Reading reference histograms for " << dirname << " from run type folder: " << foundDir->GetName());
107
108 TIter next(foundDir->GetListOfKeys());
109 TKey* hh;
110
111 while ((hh = (TKey*)next())) {
112 if (hh->IsFolder()) continue;
113 TObject* obj = hh->ReadObj(); // ReadObj -> I own it
114 if (obj->IsA()->InheritsFrom("TH1")) {
115 TH1* h = (TH1*)obj;
116 if (h->GetDimension() == 1) {
117 string histname = h->GetName();
118 m_pnode.push_back(REFNODE());
119 auto& n = m_pnode.back();
120 n.orghist_name = dirname + "/" + histname;
121 n.refhist_name = "ref/" + dirname + "/" + histname;
122 h->SetName((n.refhist_name).c_str());
123 h->SetDirectory(0);
124 n.ref_org = h; // transfer ownership!
125 n.ref_clone = nullptr;
126 n.canvas = nullptr;
127 } else {
128 delete h;
129 }
130 } else {
131 delete obj;
132 }
133 }
134 delete foundDir; // always non-zero
135 }
136 }
137 delete refdir; // always non-zero
138 }
139 }
140
141 B2INFO("DQMHistReference: insert reference to m_pnode. size: " << m_pnode.size());
142 refFile->Close();
143 delete refFile;
144}
145
147{
148 TH1::AddDirectory(false); // do not store any histograms
149
150 if (m_firstInRun) {
152 m_firstInRun = false;
153 }
154
155 char mbstr[100];
156
157 time_t now = time(0);
158 strftime(mbstr, sizeof(mbstr), "%F %T", localtime(&now));
159 B2INFO("[" << mbstr << "] before ref loop");
160
161 for (auto& it : m_pnode) {
162 if (!it.ref_org) continue; // No reference, continue
163
164 TH1* hist1 = findHistInCanvas(it.orghist_name, &(it.canvas));
165
166 TCanvas* canvas = it.canvas;
167
168 // if there is no histogram on canvas we plot the reference anyway.
169 if (!canvas) {
170 B2DEBUG(1, "No canvas found for reference histogram " << it.orghist_name);
171 continue;
172 }
173 if (!hist1) {
174 B2DEBUG(1, "Canvas is without histogram -> no display " << it.orghist_name);
175 // Display something could be confusing for shifters
176// B2DEBUG(1, "Canvas is without histogram -> displaying only reference " << it.orghist_name);
177// canvas->cd();
178// hist2->Draw();
179// canvas->Modified();
180// canvas->Update();
181 continue;
182 }
183
184 if (hist1->Integral() == 0) continue; // empty histogram -> continue
185
186 /* consider adding coloring option....
187 double data = 0;
188 if (m_color) {
189 data = hist1->KolmogorovTest(hist2, ""); // returns p value (0 bad, 1 good), N - do not compare normalized
190 }
191 */
192
193 if (abs(it.ref_org->Integral()) > 0) { // onyl if we have entries in reference
194 if (it.ref_clone) {
195 it.ref_clone->Reset();
196 it.ref_clone->Add(it.ref_org);
197 } else {
198 if (hist1->InheritsFrom("TH1C") or hist1->InheritsFrom("TH1S")) {
199 it.ref_clone = new TH1F(); // we want it a float for better scaling
200 it.ref_org->Copy(*it.ref_clone);
201 } else if (hist1->InheritsFrom("TH1I") or hist1->InheritsFrom("TH1L")) {
202 it.ref_clone = new TH1D(); // we want it a float for better scaling
203 it.ref_org->Copy(*it.ref_clone);
204 } else {
205 // keep TProfile, TH1F or TH1D
206 it.ref_clone = (TH1*)it.ref_org->Clone();
207 }
208 it.ref_clone->SetLineStyle(2);
209 it.ref_clone->SetLineColor(3);
210 it.ref_clone->SetFillColor(0);
211 it.ref_clone->SetStats(kFALSE);
212 }
213 it.ref_clone->Scale(hist1->Integral() / it.ref_clone->Integral());
214
215 // Adjust the y scale to cover the reference
216 if (it.ref_clone->GetMaximum() > hist1->GetMaximum())
217 hist1->SetMaximum(1.1 * it.ref_clone->GetMaximum());
218
219 canvas->cd();
220 it.ref_clone->Draw("hist,same");
221
222 canvas->Modified();
223 canvas->Update();
224 }
225 }
226
227 now = time(0);
228 strftime(mbstr, sizeof(mbstr), "%F %T", localtime(&now));
229 B2INFO("[" << mbstr << "] after ref loop");
230}
231
233{
234 B2DEBUG(1, "DQMHistReference: endRun called");
235}
236
237
239{
240 B2DEBUG(1, "DQMHistReference: terminate called");
241 for (auto& it : m_pnode) {
242 // clear ref histos from memory
243 if (it.ref_org) delete it.ref_org;
244 if (it.ref_clone) delete it.ref_clone;
245 }
246 m_pnode.clear();
247}
248
The base class for the histogram analysis module.
static const std::string & getRunType(void)
Get the Run Type.
TH1 * findHistInCanvas(const std::string &hname, TCanvas **canvas=nullptr)
Find histogram in corresponding canvas.
void initialize() override final
Initializer.
void loadReferenceHistos()
Reads reference histograms from input root file.
std::string m_referenceFile
Reference Histogram Root file name.
std::vector< REFNODE > m_pnode
Struct for reference histogram
bool m_firstInRun
Is first event in run.
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
void endRun() override final
This method is called if the current run ends.
void beginRun() override final
Called when entering a new run.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.
STL namespace.
Struct for refence histogram info.