Belle II Software  release-06-00-14
DQMHistAnalysisTOP.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/DQMHistAnalysisTOP.h>
10 #include <boost/format.hpp>
11 #include <daq/slc/base/StringUtil.h>
12 #include <TClass.h>
13 #include <TH2.h>
14 #include "TROOT.h"
15 
16 using namespace std;
17 using namespace Belle2;
18 using boost::format;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisTOP)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
31 {
32  //Parameter definition
33  B2DEBUG(20, "DQMHistAnalysisTOP: Constructor done.");
34 }
35 
36 
37 DQMHistAnalysisTOPModule::~DQMHistAnalysisTOPModule() { }
38 
39 void DQMHistAnalysisTOPModule::initialize()
40 {
41  gROOT->cd();
42  B2DEBUG(20, "DQMHistAnalysisTOP: initialized.");
43 
44  m_c_goodHitsMean = new TCanvas("TOP/c_good_hits_mean");
45  m_c_goodHitsRMS = new TCanvas("TOP/c_good_hits_rms");
46  m_c_badHitsMean = new TCanvas("TOP/c_bad_hits_mean");
47  m_c_badHitsRMS = new TCanvas("TOP/c_bad_hits_rms");
48 
49  //using c2_Name to avoid an overlap on default c_Name
50  for (int i = 1; i <= 16; i++) {
51  m_c_good_hits_xy_[i] = new TCanvas(Form("TOP/c2_good_hits_xy_%d", i));
52  m_c_bad_hits_xy_[i] = new TCanvas(Form("TOP/c2_bad_hits_xy_%d", i));
53  m_c_good_hits_asics_[i] = new TCanvas(Form("TOP/c2_good_hits_asics_%d", i));
54  m_c_bad_hits_asics_[i] = new TCanvas(Form("TOP/c2_bad_hits_asics_%d", i));
55  }
56 
57  m_h_goodHitsMean = new TH1F("TOP/good_hits_mean", "Mean of good hits per event", 16, 0.5, 16.5);
58  m_h_goodHitsRMS = new TH1F("TOP/good_hits_rms", "RMS of good hits per event", 16, 0.5, 16.5);
59  m_h_badHitsMean = new TH1F("TOP/bad_hits_mean", "Mean of bad hits per event", 16, 0.5, 16.5);
60  m_h_badHitsRMS = new TH1F("TOP/bad_hits_rms", "RMS of bad hits per event", 16, 0.5, 16.5);
61 
62  m_h_goodHitsMean->GetXaxis()->SetTitle("slot number");
63  m_h_goodHitsMean->GetYaxis()->SetTitle("hits per event");
64  m_h_goodHitsRMS->GetXaxis()->SetTitle("slot number");
65  m_h_goodHitsRMS->GetYaxis()->SetTitle("hits per event");
66 
67  m_h_badHitsMean->GetXaxis()->SetTitle("slot number");
68  m_h_badHitsMean->GetYaxis()->SetTitle("hits per event");
69  m_h_badHitsRMS->GetXaxis()->SetTitle("slot number");
70  m_h_badHitsRMS->GetYaxis()->SetTitle("hits per event");
71 
72  m_h_goodHitsMean->SetStats(kFALSE);
73  m_h_goodHitsRMS->SetStats(kFALSE);
74  m_h_badHitsMean->SetStats(kFALSE);
75  m_h_badHitsRMS->SetStats(kFALSE);
76 
77  m_h_goodHitsMean->SetMinimum(0);
78  m_h_goodHitsRMS->SetMinimum(0);
79  m_h_badHitsMean->SetMinimum(0);
80  m_h_badHitsRMS->SetMinimum(0);
81 
82  m_line1 = new TLine(0.5, 215, 16.5, 215);
83  m_line2 = new TLine(0.5, 235, 16.5, 235);
84  m_line1->SetLineWidth(2);
85  m_line2->SetLineWidth(2);
86  m_line1->SetLineColor(kRed);
87  m_line2->SetLineColor(kRed);
88 
89  m_text1 = new TPaveText(1, 435, 10, 500, "NB");
90  m_text1->SetFillColorAlpha(kWhite, 0);
91  m_text1->SetBorderSize(0);
92 
93  m_text2 = new TPaveText(0.2, 0.70, 0.50, 0.9, "NB");
94  m_text2->SetFillColorAlpha(kWhite, 0);
95  m_text2->SetBorderSize(0);
96 }
97 
98 
99 void DQMHistAnalysisTOPModule::beginRun()
100 {
101  //B2DEBUG(20, "DQMHistAnalysisTOP: beginRun called.");
102 }
103 
104 TH1* DQMHistAnalysisTOPModule::find_histo_in_canvas(TString histo_name)
105 {
106  StringList s = StringUtil::split(histo_name.Data(), '/');
107  std::string dirname = s[0];
108  std::string hname = s[1];
109  std::string canvas_name = dirname + "/c_" + hname;
110 
111  TIter nextckey(gROOT->GetListOfCanvases());
112 
113  auto cobj = find_canvas(canvas_name);
114  if (cobj == nullptr) return nullptr;
115 
116  TIter nextkey(((TCanvas*)cobj)->GetListOfPrimitives());
117  TObject* obj{};
118  while ((obj = dynamic_cast<TObject*>(nextkey()))) {
119  if (obj->IsA()->InheritsFrom("TH1")) {
120  if (obj->GetName() == histo_name)
121  return dynamic_cast<TH1*>(obj);
122  }
123  }
124  return nullptr;
125 }
126 
127 void DQMHistAnalysisTOPModule::event()
128 {
129  for (int i = 1; i <= 16; i++) {
130  string hname1 = str(format("TOP/good_hits_per_event%1%") % (i));;
131  string hname2 = str(format("TOP/bad_hits_per_event%1%") % (i));;
132  TH1* h1 = findHist(hname1);
133  TH1* h2 = findHist(hname2);
134  //TH1* h1 = find_histo_in_canvas(hname1);
135  //TH1* h2 = find_histo_in_canvas(hname2);
136  if (h1 != NULL) {
137  m_h_goodHitsMean->SetBinContent(i, h1->GetMean());
138  m_h_goodHitsRMS->SetBinContent(i, h1->GetRMS());
139  } else {
140  m_h_goodHitsMean->SetBinContent(i, 0);
141  m_h_goodHitsRMS->SetBinContent(i, 0);
142  }
143  if (h2 != NULL) {
144  m_h_badHitsMean->SetBinContent(i, h2->GetMean());
145  m_h_badHitsRMS->SetBinContent(i, h2->GetRMS());
146  } else {
147  m_h_badHitsMean->SetBinContent(i, 0);
148  m_h_badHitsRMS->SetBinContent(i, 0);
149  }
150  }
151 
152  TH2F* hraw = (TH2F*)findHist("TOP/window_vs_slot");
153  double exRatio(0.0);
154  double totalraw(0.0);
155  double totalbadraw(0.0);
156  if (hraw != NULL) {
157  totalraw = hraw->GetEntries();
158  for (int i = 1; i <= 16; i++) {
159  for (int j = 1; j <= 512; j++) {
160  double nhraw = hraw->GetBinContent(i, j);
161  if (j < 215 || j > 235) totalbadraw += nhraw ;
162  }
163  }
164  }
165  if (totalraw > 0) exRatio = totalbadraw * 1.0 / totalraw;
166 
167  m_text1->Clear();
168  m_text1->AddText(Form("Ratio of entries outside of red lines: %.2f %%", exRatio * 100.0));
169  if (exRatio > 0.01) {
170  m_text1->AddText(">1% bad, report to TOP experts!");
171  } else {
172  m_text1->AddText("<0.1% good, 0.1-1% recoverable.");
173  }
174 
175  //addHist("", m_h_goodHitsMean->GetName(), m_h_goodHitsMean);
176 
177  TCanvas* c1 = find_canvas("TOP/c_hitsPerEvent");
178  //TH1* h1=find_histo_in_canvas("TOP/hitsPerEvent");
179  if (c1 != NULL) {
180  c1->SetName("TOP/c_hitsPerEvent_top");
181  }
182  //if (h1!=NULL) {
183  // h1->SetName("TOP/hitsPerEvent_top");
184  //}
185 
186  TCanvas* c2 = find_canvas("TOP/c_window_vs_slot");
187  if (c2 != NULL) {
188  c2->cd();
189  if (exRatio > 0.01) c2->Pad()->SetFillColor(kRed);
190  else c2->Pad()->SetFillColor(kWhite);
191  m_line1->Draw();
192  m_line2->Draw();
193  m_text1->Draw();
194  }
195 
196  TH1F* hBoolEvtMonitor = (TH1F*)findHist("TOP/BoolEvtMonitor");
197  double badRatio(0.0);
198  int totalBadEvts(0);
199  int totalEvts(0);
200  if (hBoolEvtMonitor != NULL) {
201  totalEvts = hBoolEvtMonitor->GetEntries();
202  totalBadEvts = hBoolEvtMonitor->GetBinContent(2);
203  }
204  if (totalEvts > 0) badRatio = totalBadEvts * 1.0 / totalEvts;
205 
206  m_text2->Clear();
207  m_text2->AddText(Form("fraction of deviating hits: %.4f %%", badRatio * 100.0));
208 
209  TCanvas* c3 = find_canvas("TOP/c_BoolEvtMonitor");
210  if (c3 != NULL) {
211  c3->cd();
212  if (badRatio > 0.0001) c3->Pad()->SetFillColor(kRed);
213  else c3->Pad()->SetFillColor(kWhite);
214  m_text2->Draw();
215  }
216 
217  //obtaining the total yield for 16 2D-plots
218  double Ntotal_good_hits_xy(0.0);
219  double Ntotal_bad_hits_xy(0.0);
220  double Ntotal_good_hits_asics(0.0);
221  double Ntotal_bad_hits_asics(0.0);
222  for (int module = 1; module <= 16; module++) {
223  string hnameTmp1 = str(format("TOP/good_hits_xy_%1%") % (module));;
224  TH1* h2DTmp1 = findHist(hnameTmp1);
225  if (h2DTmp1 != NULL) Ntotal_good_hits_xy += h2DTmp1->Integral();
226 
227  string hnameTmp2 = str(format("TOP/bad_hits_xy_%1%") % (module));;
228  TH1* h2DTmp2 = findHist(hnameTmp2);
229  if (h2DTmp2 != NULL) Ntotal_bad_hits_xy += h2DTmp2->Integral();
230 
231  string hnameTmp3 = str(format("TOP/good_hits_asics_%1%") % (module));;
232  TH1* h2DTmp3 = findHist(hnameTmp3);
233  if (h2DTmp3 != NULL) Ntotal_good_hits_asics += h2DTmp3->Integral();
234 
235  string hnameTmp4 = str(format("TOP/bad_hits_asics_%1%") % (module));;
236  TH1* h2DTmp4 = findHist(hnameTmp4);
237  if (h2DTmp4 != NULL) Ntotal_bad_hits_asics += h2DTmp4->Integral();
238  }
239 
240  // reset the maximum z-axis of 16 2D plots: 3 times of average for good hits; 30 times of average for bad hits
241  for (int i = 1; i <= 16; i++) {
242  m_c_good_hits_xy_[i]->Clear();
243  m_c_good_hits_xy_[i]->cd();
244  TH2F* h2Dscale_xy = (TH2F*)findHist(Form("TOP/good_hits_xy_%d", i));
245  if (h2Dscale_xy != NULL && Ntotal_good_hits_xy > 0) {
246  h2Dscale_xy->GetZaxis()->SetRangeUser(0, Ntotal_good_hits_xy / 2500.0);
247  h2Dscale_xy->Draw();
248  }
249  m_c_good_hits_xy_[i]->Modified();
250  }
251 
252 
253  for (int i = 1; i <= 16; i++) {
254  m_c_bad_hits_xy_[i]->Clear();
255  m_c_bad_hits_xy_[i]->cd();
256  TH2F* h2Dscale_xy = (TH2F*)findHist(Form("TOP/bad_hits_xy_%d", i));
257  if (h2Dscale_xy != NULL && Ntotal_bad_hits_xy > 0) {
258  h2Dscale_xy->GetZaxis()->SetRangeUser(0, Ntotal_bad_hits_xy / 250.0);
259  h2Dscale_xy->Draw();
260  }
261  m_c_bad_hits_xy_[i]->Modified();
262  }
263 
264  for (int i = 1; i <= 16; i++) {
265  m_c_good_hits_asics_[i]->Clear();
266  m_c_good_hits_asics_[i]->cd();
267  TH2F* h2Dscale_asics = (TH2F*)findHist(Form("TOP/good_hits_asics_%d", i));
268  if (h2Dscale_asics != NULL && Ntotal_good_hits_asics > 0) {
269  h2Dscale_asics->GetZaxis()->SetRangeUser(0, Ntotal_good_hits_asics / 2500.0);
270  h2Dscale_asics->Draw();
271  }
272  m_c_good_hits_asics_[i]->Modified();
273  }
274 
275  for (int i = 1; i <= 16; i++) {
276  m_c_bad_hits_asics_[i]->Clear();
277  m_c_bad_hits_asics_[i]->cd();
278  TH2F* h2Dscale_asics = (TH2F*)findHist(Form("TOP/bad_hits_asics_%d", i));
279  if (h2Dscale_asics != NULL && Ntotal_bad_hits_asics > 0) {
280  h2Dscale_asics->GetZaxis()->SetRangeUser(0, Ntotal_bad_hits_asics / 250.0);
281  h2Dscale_asics->Draw();
282  }
283  m_c_bad_hits_asics_[i]->Modified();
284  }
285 
286  m_c_goodHitsMean->Clear();
287  m_c_goodHitsMean->cd();
288  m_h_goodHitsMean->Draw();
289  m_c_goodHitsMean->Modified();
290 
291  m_c_goodHitsRMS->Clear();
292  m_c_goodHitsRMS->cd();
293  m_h_goodHitsRMS->Draw();
294  m_c_goodHitsRMS->Modified();
295 
296  m_c_badHitsMean->Clear();
297  m_c_badHitsMean->cd();
298  m_h_badHitsMean->Draw();
299  m_c_badHitsMean->Modified();
300 
301  m_c_badHitsRMS->Clear();
302  m_c_badHitsRMS->cd();
303  m_h_badHitsRMS->Draw();
304  m_c_badHitsRMS->Modified();
305 }
306 
307 void DQMHistAnalysisTOPModule::endRun()
308 {
309  B2DEBUG(20, "DQMHistAnalysisTOP : endRun called");
310 }
311 
312 
313 void DQMHistAnalysisTOPModule::terminate()
314 {
315  B2DEBUG(20, "terminate called");
316 }
317 
The base class for the histogram analysis module.
Class for TOP histogram analysis.
#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.