Belle II Software  release-05-01-25
DQMHistAnalysisTOP.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2018 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Boqun Wang *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <dqm/analysis/modules/DQMHistAnalysisTOP.h>
12 #include <boost/format.hpp>
13 #include <daq/slc/base/StringUtil.h>
14 #include <TClass.h>
15 #include <TH2.h>
16 #include "TROOT.h"
17 
18 using namespace std;
19 using namespace Belle2;
20 using boost::format;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(DQMHistAnalysisTOP)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
33 {
34  //Parameter definition
35  B2DEBUG(20, "DQMHistAnalysisTOP: Constructor done.");
36 }
37 
38 
39 DQMHistAnalysisTOPModule::~DQMHistAnalysisTOPModule() { }
40 
41 void DQMHistAnalysisTOPModule::initialize()
42 {
43  gROOT->cd();
44  B2DEBUG(20, "DQMHistAnalysisTOP: initialized.");
45 
46  m_c_goodHitsMean = new TCanvas("TOP/c_good_hits_mean");
47  m_c_goodHitsRMS = new TCanvas("TOP/c_good_hits_rms");
48  m_c_badHitsMean = new TCanvas("TOP/c_bad_hits_mean");
49  m_c_badHitsRMS = new TCanvas("TOP/c_bad_hits_rms");
50 
51  m_h_goodHitsMean = new TH1F("TOP/good_hits_mean", "Mean of good hits per event", 16, 0.5, 16.5);
52  m_h_goodHitsRMS = new TH1F("TOP/good_hits_rms", "RMS of good hits per event", 16, 0.5, 16.5);
53  m_h_badHitsMean = new TH1F("TOP/bad_hits_mean", "Mean of bad hits per event", 16, 0.5, 16.5);
54  m_h_badHitsRMS = new TH1F("TOP/bad_hits_rms", "RMS of bad hits per event", 16, 0.5, 16.5);
55 
56  m_h_goodHitsMean->GetXaxis()->SetTitle("slot number");
57  m_h_goodHitsMean->GetYaxis()->SetTitle("hits per event");
58  m_h_goodHitsRMS->GetXaxis()->SetTitle("slot number");
59  m_h_goodHitsRMS->GetYaxis()->SetTitle("hits per event");
60 
61  m_h_badHitsMean->GetXaxis()->SetTitle("slot number");
62  m_h_badHitsMean->GetYaxis()->SetTitle("hits per event");
63  m_h_badHitsRMS->GetXaxis()->SetTitle("slot number");
64  m_h_badHitsRMS->GetYaxis()->SetTitle("hits per event");
65 
66  m_h_goodHitsMean->SetStats(kFALSE);
67  m_h_goodHitsRMS->SetStats(kFALSE);
68  m_h_badHitsMean->SetStats(kFALSE);
69  m_h_badHitsRMS->SetStats(kFALSE);
70 
71  m_h_goodHitsMean->SetMinimum(0);
72  m_h_goodHitsRMS->SetMinimum(0);
73  m_h_badHitsMean->SetMinimum(0);
74  m_h_badHitsRMS->SetMinimum(0);
75 
76  m_line1 = new TLine(0.5, 215, 16.5, 215);
77  m_line2 = new TLine(0.5, 235, 16.5, 235);
78  m_line1->SetLineWidth(2);
79  m_line2->SetLineWidth(2);
80  m_line1->SetLineColor(kRed);
81  m_line2->SetLineColor(kRed);
82 
83  m_text1 = new TPaveText(1, 435, 10, 500, "NB");
84  m_text1->SetFillColorAlpha(kWhite, 0);
85  m_text1->SetBorderSize(0);
86 
87  m_text2 = new TPaveText(0.2, 0.70, 0.50, 0.9, "NB");
88  m_text2->SetFillColorAlpha(kWhite, 0);
89  m_text2->SetBorderSize(0);
90 }
91 
92 
93 void DQMHistAnalysisTOPModule::beginRun()
94 {
95  //B2DEBUG(20, "DQMHistAnalysisTOP: beginRun called.");
96 }
97 
98 TCanvas* DQMHistAnalysisTOPModule::find_canvas(TString canvas_name)
99 {
100  TIter nextckey(gROOT->GetListOfCanvases());
101  TObject* cobj = NULL;
102 
103  while ((cobj = (TObject*)nextckey())) {
104  if (cobj->IsA()->InheritsFrom("TCanvas")) {
105  if (cobj->GetName() == canvas_name)
106  break;
107  }
108  }
109  return (TCanvas*)cobj;
110 }
111 
112 TH1* DQMHistAnalysisTOPModule::find_histo_in_canvas(TString histo_name)
113 {
114  StringList s = StringUtil::split(histo_name.Data(), '/');
115  std::string dirname = s[0];
116  std::string hname = s[1];
117  std::string canvas_name = dirname + "/c_" + hname;
118 
119  TIter nextckey(gROOT->GetListOfCanvases());
120  TObject* cobj = NULL;
121 
122  while ((cobj = (TObject*)nextckey())) {
123  if (cobj->IsA()->InheritsFrom("TCanvas")) {
124  if (cobj->GetName() == canvas_name)
125  break;
126  }
127  }
128  if (cobj == NULL) return NULL;
129 
130  TIter nextkey(((TCanvas*)cobj)->GetListOfPrimitives());
131  TObject* obj = NULL;
132 
133  while ((obj = (TObject*)nextkey())) {
134  if (obj->IsA()->InheritsFrom("TH1")) {
135  if (obj->GetName() == histo_name)
136  return (TH1*)obj;
137  }
138  }
139  return NULL;
140 }
141 
142 void DQMHistAnalysisTOPModule::event()
143 {
144  for (int i = 1; i <= 16; i++) {
145  string hname1 = str(format("TOP/good_hits_per_event%1%") % (i));;
146  string hname2 = str(format("TOP/bad_hits_per_event%1%") % (i));;
147  TH1* h1 = findHist(hname1);
148  TH1* h2 = findHist(hname2);
149  //TH1* h1 = find_histo_in_canvas(hname1);
150  //TH1* h2 = find_histo_in_canvas(hname2);
151  if (h1 != NULL) {
152  m_h_goodHitsMean->SetBinContent(i, h1->GetMean());
153  m_h_goodHitsRMS->SetBinContent(i, h1->GetRMS());
154  } else {
155  m_h_goodHitsMean->SetBinContent(i, 0);
156  m_h_goodHitsRMS->SetBinContent(i, 0);
157  }
158  if (h2 != NULL) {
159  m_h_badHitsMean->SetBinContent(i, h2->GetMean());
160  m_h_badHitsRMS->SetBinContent(i, h2->GetRMS());
161  } else {
162  m_h_badHitsMean->SetBinContent(i, 0);
163  m_h_badHitsRMS->SetBinContent(i, 0);
164  }
165  }
166 
167  TH2F* hraw = (TH2F*)findHist("TOP/window_vs_slot");
168  double exRatio(0.0);
169  double totalraw(0.0);
170  double totalbadraw(0.0);
171  if (hraw != NULL) {
172  totalraw = hraw->GetEntries();
173  for (int i = 1; i <= 16; i++) {
174  for (int j = 1; j <= 512; j++) {
175  double nhraw = hraw->GetBinContent(i, j);
176  if (j < 215 || j > 235) totalbadraw += nhraw ;
177  }
178  }
179  }
180  if (totalraw > 0) exRatio = totalbadraw * 1.0 / totalraw;
181 
182  m_text1->Clear();
183  m_text1->AddText(Form("Ratio of entries outside of red lines: %.2f %%", exRatio * 100.0));
184  if (exRatio > 0.01) {
185  m_text1->AddText(">1% bad, report to TOP experts!");
186  } else {
187  m_text1->AddText("<0.1% good, 0.1-1% recoverable.");
188  }
189 
190  //addHist("", m_h_goodHitsMean->GetName(), m_h_goodHitsMean);
191 
192  TCanvas* c1 = find_canvas("TOP/c_hitsPerEvent");
193  //TH1* h1=find_histo_in_canvas("TOP/hitsPerEvent");
194  if (c1 != NULL) {
195  c1->SetName("TOP/c_hitsPerEvent_top");
196  }
197  //if (h1!=NULL) {
198  // h1->SetName("TOP/hitsPerEvent_top");
199  //}
200 
201  TCanvas* c2 = find_canvas("TOP/c_window_vs_slot");
202  if (c2 != NULL) {
203  c2->cd();
204  if (exRatio > 0.01) c2->Pad()->SetFillColor(kRed);
205  else c2->Pad()->SetFillColor(kWhite);
206  m_line1->Draw();
207  m_line2->Draw();
208  m_text1->Draw();
209  }
210 
211  TH1F* hBoolEvtMonitor = (TH1F*)findHist("TOP/BoolEvtMonitor");
212  double badRatio(0.0);
213  int totalBadEvts(0);
214  int totalEvts(0);
215  if (hBoolEvtMonitor != NULL) {
216  totalEvts = hBoolEvtMonitor->GetEntries();
217  totalBadEvts = hBoolEvtMonitor->GetBinContent(2);
218  }
219  if (totalEvts > 0) badRatio = totalBadEvts * 1.0 / totalEvts;
220 
221  m_text2->Clear();
222  m_text2->AddText(Form("fraction of deviating hits: %.4f %%", badRatio * 100.0));
223 
224  TCanvas* c3 = find_canvas("TOP/c_BoolEvtMonitor");
225  if (c3 != NULL) {
226  c3->cd();
227  if (badRatio > 0.0001) c3->Pad()->SetFillColor(kRed);
228  else c3->Pad()->SetFillColor(kWhite);
229  m_text2->Draw();
230  }
231 
232  m_c_goodHitsMean->Clear();
233  m_c_goodHitsMean->cd();
234  m_h_goodHitsMean->Draw();
235  m_c_goodHitsMean->Modified();
236 
237  m_c_goodHitsRMS->Clear();
238  m_c_goodHitsRMS->cd();
239  m_h_goodHitsRMS->Draw();
240  m_c_goodHitsRMS->Modified();
241 
242  m_c_badHitsMean->Clear();
243  m_c_badHitsMean->cd();
244  m_h_badHitsMean->Draw();
245  m_c_badHitsMean->Modified();
246 
247  m_c_badHitsRMS->Clear();
248  m_c_badHitsRMS->cd();
249  m_h_badHitsRMS->Draw();
250  m_c_badHitsRMS->Modified();
251 }
252 
253 void DQMHistAnalysisTOPModule::endRun()
254 {
255  B2DEBUG(20, "DQMHistAnalysisTOP : endRun called");
256 }
257 
258 
259 void DQMHistAnalysisTOPModule::terminate()
260 {
261  B2DEBUG(20, "terminate called");
262 }
263 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisTOPModule
Class for TOP histogram analysis.
Definition: DQMHistAnalysisTOP.h:38
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27