Belle II Software  release-08-01-10
DQMHistComparitor.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 : DQMHistComparitor.cc
10 // Description : DQM Histogram analysis module, compares histo with references
11 //-
12 
13 
14 #include <framework/core/ModuleParam.templateDetails.h>
15 #include <dqm/analysis/modules/DQMHistComparitor.h>
16 #include <TROOT.h>
17 #include <TStyle.h>
18 #include <TClass.h>
19 #include <TDirectory.h>
20 #include <TH1F.h>
21 #include <TH2F.h>
22 #include <TKey.h>
23 
24 using namespace std;
25 using namespace Belle2;
26 
27 //-----------------------------------------------------------------
28 // Register the Module
29 //-----------------------------------------------------------------
30 REG_MODULE(DQMHistComparitor);
31 
32 //-----------------------------------------------------------------
33 // Implementation
34 //-----------------------------------------------------------------
35 
36 DQMHistComparitorModule::DQMHistComparitorModule()
38 {
39  //Parameter definition
40  addParam("HistoList", m_histlist, "['histname', 'refhistname', 'canvas', 'minentry', 'warnlvl', 'errlvl', 'pvname'],[...],...");
41  addParam("RefHistoFile", m_refFileName, "Reference histrogram file name", std::string("refHisto.root"));
42  addParam("ColorAlert", m_color, "Whether to show the color alert", true);
43  B2DEBUG(1, "DQMHistComparitor: Constructor done.");
44 }
45 
47 {
48 #ifdef _BELLE2_EPICS
49  if (ca_current_context()) ca_context_destroy();
50 #endif
51 }
52 
53 TH1* DQMHistComparitorModule::GetHisto(TString histoname)
54 {
55  TH1* hh1;
56  gROOT->cd();
57  hh1 = findHist(histoname.Data());
58  if (hh1 == NULL) {
59  B2DEBUG(20, "findHist failed " << histoname << " not in memfile");
60 
61  // first search reference root file ... if ther is one
62  if (m_refFile && m_refFile->IsOpen()) {
63  TDirectory* d = m_refFile;
64  TString myl = histoname;
65  TString tok;
66  Ssiz_t from = 0;
67  B2DEBUG(20, myl);
68  while (myl.Tokenize(tok, from, "/")) {
69  TString dummy;
70  Ssiz_t f;
71  f = from;
72  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
73  auto e = d->GetDirectory(tok);
74  if (e) {
75  B2DEBUG(20, "Cd Dir " << tok << " from " << d->GetPath());
76  d = e;
77  } else {
78  B2DEBUG(20, "cd failed " << tok << " from " << d->GetPath());
79  }
80  } else {
81  break;
82  }
83  }
84  TObject* obj = d->FindObject(tok);
85  if (obj != NULL) {
86  if (obj->IsA()->InheritsFrom("TH1")) {
87  B2DEBUG(20, "Histo " << histoname << " found in ref file");
88  hh1 = (TH1*)obj;
89  } else {
90  B2DEBUG(20, "Histo " << histoname << " found in ref file but wrong type");
91  }
92  } else {
93  // seems find will only find objects, not keys, thus get the object on first access
94  TIter next(d->GetListOfKeys());
95  TKey* key;
96  while ((key = (TKey*)next())) {
97  TObject* obj2 = key->ReadObj() ;
98  if (obj2->InheritsFrom("TH1")) {
99  if (obj2->GetName() == tok) {
100  hh1 = (TH1*)obj2;
101  B2DEBUG(20, "Histo " << histoname << " found as key -> readobj");
102  break;
103  }
104  }
105  }
106  if (hh1 == NULL) B2DEBUG(20, "Histo " << histoname << " NOT found in ref file " << tok);
107  }
108  }
109 
110  if (hh1 == NULL) {
111  B2DEBUG(20, "Histo " << histoname << " not in memfile or ref file");
112 
113  TDirectory* d = gROOT;
114  TString myl = histoname;
115  TString tok;
116  Ssiz_t from = 0;
117  while (myl.Tokenize(tok, from, "/")) {
118  TString dummy;
119  Ssiz_t f;
120  f = from;
121  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
122  auto e = d->GetDirectory(tok);
123  if (e) {
124  B2DEBUG(20, "Cd Dir " << tok);
125  d = e;
126  } else B2DEBUG(20, "cd failed " << tok);
127  d->cd();
128  } else {
129  break;
130  }
131  }
132  TObject* obj = d->FindObject(tok);
133  if (obj != NULL) {
134  if (obj->IsA()->InheritsFrom("TH1")) {
135  B2DEBUG(20, "Histo " << histoname << " found in mem");
136  hh1 = (TH1*)obj;
137  }
138  } else {
139  B2DEBUG(20, "Histo " << histoname << " NOT found in mem");
140  }
141  }
142  }
143 
144  if (hh1 == NULL) {
145  B2DEBUG(20, "Histo " << histoname << " not found");
146  }
147 
148  return hh1;
149 }
150 
152 {
153  m_refFile = NULL;
154  if (m_refFileName != "") {
155  m_refFile = new TFile(m_refFileName.data());
156  }
157 
158  gStyle->SetOptStat(0);
159  gStyle->SetStatStyle(1);
160  gStyle->SetOptDate(22);// Date and Time in Bottom Right, does no work
161 
162 #ifdef _BELLE2_EPICS
163  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
164 #endif
165  for (auto& it : m_histlist) {
166  if (it.size() != 7) {
167  B2WARNING("Histolist with wrong nr of parameters (" << it.size() << "), hist1=" << it.at(0));
168  continue;
169  }
170 
171  // Do not check for histogram existance here, as it might not be
172  // in memfile when analysis task is started
173  TString a = it.at(2).c_str();
174 
175  auto n = new CMPNODE;
176  n->histo1 = it.at(0);
177  n->histo2 = it.at(1);
178  TCanvas* c = new TCanvas(a);
179  n->canvas = c;
180 
181  n->min_entries = atoi(it.at(3).c_str());
182  n->warning = atof(it.at(4).c_str());
183  n->error = atof(it.at(5).c_str());
184  n->epicsflag = false;
185 
186 #ifdef _BELLE2_EPICS
187  if (it.at(6) != "") {
188  SEVCHK(ca_create_channel(it.at(6).c_str(), NULL, NULL, 10, &n->mychid), "ca_create_channel failure");
189  n->epicsflag = true;
190  }
191 #endif
192  m_pnode.push_back(n);
193  }
194 #ifdef _BELLE2_EPICS
195  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
196 #endif
197 
198  B2DEBUG(20, "DQMHistComparitor: initialized.");
199 }
200 
201 
203 {
204  B2DEBUG(20, "DQMHistComparitor: beginRun called.");
205 }
206 
208 {
209  for (auto& it : m_pnode) {
210 
211  TH1* hist1, *hist2;
212 
213  B2DEBUG(20, "== Search for " << it->histo1 << " with ref " << it->histo2 << "==");
214  hist1 = findHistInCanvas(it->histo1.Data());
215  if (!hist1) B2DEBUG(20, "NOT Found " << it->histo1);
216  hist2 = GetHisto(it->histo2);
217  if (!hist2) B2DEBUG(20, "NOT Found " << it->histo2);
218  // Dont compare if any of hists is missing ... TODO We should clean CANVAS, raise some error if we cannot compare
219  if (!hist1) continue;
220  if (!hist2) continue;
221 
222  B2DEBUG(20, "Compare " << it->histo1 << " with ref " << it->histo2);
223  // if compare normalized ... does not work!
224  // hist2->Scale(hist1->GetEntries()/hist2->GetEntries());
225 
226  double data = hist1->KolmogorovTest(hist2, ""); // returns p value (0 bad, 1 good), N - do not compare normalized
227  // data = hist1->Chi2Test(hist2);// return p value (0 bad, 1 good), ignores normalization
228  // data= BinByBinTest(hits1,hist2);// user function (like Peters test)
229  // printf(" %.2f %.2f %.2f\n",(float)data,it->warning,it->error);
230 #ifdef _BELLE2_EPICS
231  if (it->epicsflag) SEVCHK(ca_put(DBR_DOUBLE, it->mychid, (void*)&data), "ca_set failure");
232 #endif
233  it->canvas->cd();
234  hist2->SetLineStyle(3);// 2 or 3
235  hist2->SetLineColor(3);
236 
237  TIter nextkey(it->canvas->GetListOfPrimitives());
238  TObject* obj = NULL;
239  while ((obj = (TObject*)nextkey())) {
240  if (obj->IsA()->InheritsFrom("TH1")) {
241  if (string(obj->GetName()) == string(hist2->GetName())) {
242  delete obj;
243  }
244  }
245  }
246 
247  // if draw normalized
248  TH1* h;
249  if (1) {
250  h = (TH1*)hist2->Clone(); // Anoying ... Maybe an memory leak? TODO
251  if (abs(hist2->Integral()) > 0)
252  h->Scale(hist1->Integral() / hist2->Integral());
253  } else {
254  hist2->Draw("hist");
255  }
256 
257  h->SetFillColor(0);
258  h->SetStats(kFALSE);
259  hist1->SetFillColor(0);
260  if (h->GetMaximum() > hist1->GetMaximum())
261  hist1->SetMaximum(1.1 * h->GetMaximum());
262  hist1->Draw("hist");
263  h->Draw("hist,same");
264 
265  it->canvas->Pad()->SetFrameFillColor(10);
266  if (m_color) {
267  if (hist1->GetEntries() < it->min_entries) {
268  // not enough Entries
269  it->canvas->Pad()->SetFillColor(6);// Magenta
270  } else {
271  if (data < it->error) {
272  it->canvas->Pad()->SetFillColor(2);// Red
273  } else if (data < it->warning) {
274  it->canvas->Pad()->SetFillColor(5);// Yellow
275  } else {
276  it->canvas->Pad()->SetFillColor(0);// White
277  }
278  }
279  } else {
280  it->canvas->Pad()->SetFillColor(0);// White
281  }
282  it->canvas->Modified();
283  it->canvas->Update();
284  }
285 #ifdef _BELLE2_EPICS
286  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
287 #endif
288 }
289 
291 {
292  B2DEBUG(20, "DQMHistComparitor: endRun called");
293 }
294 
295 
297 {
298  B2DEBUG(20, "DQMHistComparitor: terminate called");
299  if (m_refFile) delete m_refFile;
300 }
The base class for the histogram analysis module.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
TH1 * findHistInCanvas(const std::string &hname)
Find histogram in corresponding canvas.
std::string m_refFileName
Reference Histogram Root file name.
bool m_color
Whether to use the color code for warnings and errors.
void initialize() override final
Initializer.
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.
TH1 * GetHisto(TString histoname)
Get histogram by its name.
void endRun() override final
This method is called if the current run ends.
std::vector< std::vector< std::string > > m_histlist
Parameter list for histograms.
void beginRun() override final
Called when entering a new run.
std::vector< CMPNODE * > m_pnode
Struct for extracted parameters + EPICS PV.
TFile * m_refFile
The pointer to the reference file.
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.
The struct for reference histogram comparison.