Belle II Software  release-05-01-25
DQMHistAnalysisECL.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2018 - Belle II Collaboration *
4  * *
5  * ECL Data Quality Monitor (Analysis part) *
6  * *
7  * This module provides analysis of ECL DQM histograms *
8  * @ 2D Histos comparitor (color alert) *
9  * @ Histo title renew (adc_flag) *
10  * @ 2D histo analysis (time_crate_%1%_Thr1GeV) *
11  * *
12  * Author: The Belle II Collaboration *
13  * Contributors: Dmitry Matvienko (d.v.matvienko@inp.nsk.su) *
14  * *
15  * This software is provided "as is" without any warranty. *
16  **************************************************************************/
17 
18 //THIS MODULE
19 #include <dqm/analysis/modules/DQMHistAnalysisECL.h>
20 
21 #include <boost/format.hpp>
22 #include <cmath>
23 #include <TROOT.h>
24 #include <TStyle.h>
25 #include <TColor.h>
26 #include <sstream>
27 
28 using namespace Belle2;
29 
30 REG_MODULE(DQMHistAnalysisECL)
31 
34 {
35  B2DEBUG(20, "DQMHistAnalysisECL: Constructor done.");
36 
37  m_WaveformOption = {"psd", "logic", "rand", "dphy", "other"};
38 
39  addParam("HitMapThresholds", m_HitMapThresholds, "Thresholds to display hit map, MeV", std::vector<double> {0, 5, 10, 50});
40  addParam("WaveformOption", m_WaveformOption, "Option (all,psd,logic,rand,dphy) to display waveform flow",
41  m_WaveformOption);
42  addParam("CrateTimeOffsetsMax", m_CrateTimeOffsetsMax, "Maximum boundary for crate time offsets", 20.);
43  addParam("LogicTestMax", m_LogicTestMax, " Maximum of fails for logic test", 50);
44 }
45 
46 
48 
50 {
51  gROOT->cd();
52  B2DEBUG(20, "DQMHistAnalysisECL: initialized.");
53 
54  //new canvases for existing histograms
55  c_quality_analysis = new TCanvas("ECL/c_quality_analysis");
56  c_quality_other_analysis = new TCanvas("ECL/c_quality_other_analysis");
57  c_bad_quality_analysis = new TCanvas("ECL/c_bad_quality_analysis");
58  c_trigtag1_analysis = new TCanvas("ECL/c_trigtag1_analysis");
59  c_trigtag2_analysis = new TCanvas("ECL/c_trigtag2_analysis");
60  c_adc_hits_analysis = new TCanvas("ECL/c_adc_hits_analysis");
61  c_ampfail_quality_analysis = new TCanvas("ECL/c_ampfail_quality_analysis");
62  c_timefail_quality_analysis = new TCanvas("ECL/c_timefail_quality_analysis");
63  c_quality_fit_data_analysis = new TCanvas("ECL/c_quality_fit_data_analysis");
64 
65  for (const auto& id : m_HitMapThresholds) {
66  std::string canvas_name = str(boost::format("ECL/c_cid_Thr%1%MeV_analysis") % id);
67  TCanvas* canvas = new TCanvas(canvas_name.c_str());
68  c_cid_analysis.push_back(canvas);
69  }
70 
71  for (const auto& id : m_WaveformOption) {
72  if (id != "other") {
73  std::string canvas_name = str(boost::format("ECL/c_wf_cid_%1%_analysis") % id);
74  TCanvas* canvas = new TCanvas(canvas_name.c_str());
75  c_wf_analysis.push_back(canvas);
76  }
77  }
78 
79  //Boundaries for 'trigtag2_trigid' histogram
80  m_lower_boundary_trigtag2 = new TLine(1, 0, 53, 0);
81  m_lower_boundary_trigtag2->SetLineWidth(3);
82  m_lower_boundary_trigtag2->SetLineColor(kBlue);
83 
84  m_upper_boundary_trigtag2 = new TLine(1, 1, 53, 1);
85  m_upper_boundary_trigtag2->SetLineWidth(3);
86  m_upper_boundary_trigtag2->SetLineColor(kBlue);
87  //Boundaries for 'crate_time_offset' plot
88  m_lower_boundary_time_offsets = new TLine(0, -20, 52, -20);
89  m_lower_boundary_time_offsets->SetLineWidth(3);
90  m_lower_boundary_time_offsets->SetLineColor(kBlue);
91 
92  m_upper_boundary_time_offsets = new TLine(0, 20, 52, 20);
93  m_upper_boundary_time_offsets->SetLineWidth(3);
94  m_upper_boundary_time_offsets->SetLineColor(kBlue);
95 
96  //Summary crate_time_offsets plot
97  c_crate_time_offsets = new TCanvas("ECL/c_crate_time_offsets");
98  h_crate_time_offsets = new TGraphErrors();
99  h_crate_time_offsets->SetName("t_off");
100  h_crate_time_offsets->SetTitle("Crate time offset (E > 1 GeV); Crate ID (same as ECLCollector ID); Time offset [ns]");
101  h_crate_time_offsets->SetMarkerColor(kBlue);
102  h_crate_time_offsets->SetMarkerStyle(20);
103 
104  //New DQM summary for logic test in CR shifter panel
105  c_logic_summary = new TCanvas("ECL/c_logic_summary");
106  h_logic_summary = new TH2F("logic_summary", "FPGA <-> C++ fitter inconsistencies count", 52, 1, 53, 12, 1, 13);
107  h_logic_summary->SetTitle("FPGA <-> C++ fitter inconsistencies count; ECLCollector ID (same as Crate ID); Shaper ID inside the crate");
108  h_logic_summary->SetCanExtend(TH1::kAllAxes);
109  h_logic_summary->SetStats(0);
110  for (unsigned short i = 0; i < 52; i++) h_logic_summary->GetXaxis()->SetBinLabel(i + 1, std::to_string(i + 1).c_str());
111  for (unsigned short i = 0; i < 12; i++) h_logic_summary->GetYaxis()->SetBinLabel(i + 1, std::to_string(i + 1).c_str());
112  h_logic_summary->LabelsOption("v");
113  h_logic_summary->SetTickLength(0, "xy");
114 }
115 
117 {
118  B2DEBUG(20, "DQMHistAnalysisECL: beginRun called.");
119 }
120 
121 void DQMHistAnalysisECLModule::normalize(TCanvas* c, const std::string& h_name, const Double_t& weight)
122 {
123  c->Clear();
124  c->cd();
125  TH1* h = findHist(h_name);
126  if (h != NULL) {
127  for (unsigned short i = 0; i < h->GetNbinsX(); i++) {
128  Double_t entries = h->GetBinContent(i + 1);
129  Double_t error = h->GetBinError(i + 1);
130  if (weight) {
131  h->SetBinContent(i + 1, entries / weight);
132  h->SetBinError(i + 1, error / weight);
133  }
134  }
135  h->Draw("HIST");
136  }
137  c->Draw();
138  c->Modified();
139  c->Update();
140 }
141 
143 {
144  B2DEBUG(20, "DQMHistAnalysisECL: event called");
145 
146  //gStyle requirements
147  gStyle->SetPalette(1);
148 
149  //1D histos
150  //quality
151  c_quality_analysis->Clear();
152  c_quality_analysis->cd();
153  c_quality_analysis->SetLogy();
154  TH1* h_quality = findHist("ECL/quality");
155  if (h_quality != NULL) {
156  h_quality->SetMinimum(0.1);
157  h_quality->Draw();
158  }
159  c_quality_analysis->Draw();
160  c_quality_analysis->Modified();
161  c_quality_analysis->Update();
162 
163  //quality_other
164  c_quality_other_analysis->Clear();
166  c_quality_other_analysis->SetLogy();
167  TH1* h_quality_other = findHist("ECL/quality_other");
168  if (h_quality_other != NULL) {
169  h_quality_other->SetMinimum(0.1);
170  h_quality_other->Draw();
171  }
172  c_quality_other_analysis->Draw();
173  c_quality_other_analysis->Modified();
174  c_quality_other_analysis->Update();
175 
176  //bad_quality ,cid_Thr%1%MeV, wf_cid_%1%, wf_sh_%1%, wf_cr_%1%
177  TH1* h_evtot = findHist("ECL/event");
178  if (h_evtot != NULL) {
179  Double_t events = h_evtot->GetBinContent(1);
180  normalize(c_bad_quality_analysis, "ECL/bad_quality", events);
181  for (std::size_t i = 0; i < m_HitMapThresholds.size(); ++i)
183  str(boost::format("ECL/cid_Thr%1%MeV") % m_HitMapThresholds[i]), events);
184  }
185  for (std::size_t i = 0; i < m_WaveformOption.size(); ++i) {
186  auto val = m_WaveformOption[i];
187  if (val != "psd" && val != "other") {
188  TH1* h_evtot_norm = findHist(str(boost::format("ECL/event_%1%") % val));
189  if (h_evtot_norm != NULL) {
190  Double_t events = h_evtot_norm->GetBinContent(1);
192  str(boost::format("ECL/wf_cid_%1%") % val), events);
193  }
194  } else if (val == "psd") {
195  c_wf_analysis[i]->Clear();
196  c_wf_analysis[i]->cd();
197  TH1* h_psd = findHist(str(boost::format("ECL/wf_cid_%1%") % val));
198  TH1* h_psd_norm = findHist(str(boost::format("ECL/%1%_cid") % val));
199  if (h_psd != NULL && h_psd_norm != NULL) {
200  h_psd->Divide(h_psd, h_psd_norm);
201  h_psd->Draw("HIST");
202  }
203  c_wf_analysis[i]->Draw();
204  c_wf_analysis[i]->Modified();
205  c_wf_analysis[i]->Update();
206  }
207  } //m_WaveformOption
208 
209  //trigtag1
210  c_trigtag1_analysis->Clear();
211  c_trigtag1_analysis->cd();
212  c_trigtag1_analysis->Pad()->SetFrameFillColor(10);
213  c_trigtag1_analysis->Pad()->SetFillColor(kWhite);
214  c_trigtag1_analysis->SetLogy();
215  TH1* h_trigtag1 = findHist("ECL/trigtag1");
216  if (h_trigtag1 != NULL) {
217  h_trigtag1->SetMinimum(0.1);
218  h_trigtag1->Draw();
219  if (h_trigtag1->GetBinContent(2)) c_trigtag1_analysis->Pad()->SetFillColor(kRed);
220  }
221  c_trigtag1_analysis->Draw();
222  c_trigtag1_analysis->Modified();
223  c_trigtag1_analysis->Update();
224 
225  //adc_hits
226  c_adc_hits_analysis->Clear();
227  c_adc_hits_analysis->cd();
228  c_adc_hits_analysis->SetLogy();
229  TH1* h_adc_hits = findHist("ECL/adc_hits");
230  if (h_adc_hits != NULL) {
231  h_adc_hits->SetMinimum(0.1);
232  h_adc_hits->Draw();
233  }
234  c_adc_hits_analysis->Draw();
235  c_adc_hits_analysis->Modified();
236  c_adc_hits_analysis->Update();
237 
238  //ampfail_quality
241  c_ampfail_quality_analysis->Pad()->SetFrameFillColor(10);
242  c_ampfail_quality_analysis->Pad()->SetFillColor(kWhite);
243  c_ampfail_quality_analysis->SetLogy();
244  TH1* h_ampfail_quality = findHist("ECL/ampfail_quality");
245  if (h_ampfail_quality != NULL) {
246  h_ampfail_quality->SetMinimum(0.1);
247  h_ampfail_quality->Draw();
248  for (unsigned short i = 1; i < 5; i++) {
249  if (h_ampfail_quality->GetBinContent(i + 1)) {
250  c_ampfail_quality_analysis->Pad()->SetFillColor(kRed);
251  break;
252  }
253  }
254  }
256  c_ampfail_quality_analysis->Modified();
257  c_ampfail_quality_analysis->Update();
258 
259  //timefail_quality
262  c_timefail_quality_analysis->Pad()->SetFrameFillColor(10);
263  c_timefail_quality_analysis->Pad()->SetFillColor(kWhite);
264  c_timefail_quality_analysis->SetLogy();
265  TH1* h_timefail_quality = findHist("ECL/timefail_quality");
266  if (h_timefail_quality != NULL) {
267  h_timefail_quality->SetMinimum(0.1);
268  h_timefail_quality->Draw();
269  for (unsigned short i = 1; i < 5; i++) {
270  if (h_timefail_quality->GetBinContent(i + 1)) {
271  c_timefail_quality_analysis->Pad()->SetFillColor(kRed);
272  break;
273  }
274  }
275  }
277  c_timefail_quality_analysis->Modified();
278  c_timefail_quality_analysis->Update();
279 
280  //2D histos
281  //trigtag2_trigid
282  c_trigtag2_analysis->Clear();
283  c_trigtag2_analysis->cd();
284  c_trigtag2_analysis->Pad()->SetFrameFillColor(10);
285  c_trigtag2_analysis->Pad()->SetFillColor(kWhite);
286  TH1* h_trigtag2_trigid = findHist("ECL/trigtag2_trigid");
287  if (h_trigtag2_trigid != NULL) {
288  h_trigtag2_trigid->Draw();
289  for (unsigned short i = 0; i < 52; i++) {
290  if (h_trigtag2_trigid->GetBinContent(h_trigtag2_trigid->GetBin(i + 1, 3))) {
291  c_trigtag2_analysis->Pad()->SetFillColor(kRed);
292  break;
293  }
294  }
295  }
298 
299  c_trigtag2_analysis->Draw();
300  c_trigtag2_analysis->Modified();
301  c_trigtag2_analysis->Update();
302 
303  //quality_fit_data
306  c_quality_fit_data_analysis->Pad()->SetFrameFillColor(10);
307  c_quality_fit_data_analysis->Pad()->SetFillColor(kWhite);
308  TH1* h_quality_fit_data = findHist("ECL/quality_fit_data");
309  if (h_quality_fit_data != NULL) {
310  h_quality_fit_data->Draw();
311  for (unsigned short i = 0; i < 4; i++) {
312  for (unsigned short j = 0; j < 4; j++) {
313  if (h_quality_fit_data->GetBinContent(h_quality_fit_data->GetBin(i + 1, j + 1)) > 0) {
314  c_quality_fit_data_analysis->Pad()->SetFillColor(kRed);
315  break;
316  }
317  }
318  }
319  }
321  c_quality_fit_data_analysis->Modified();
322  c_quality_fit_data_analysis->Update();
323 
324  //_time_crate_%1%_Thr1GeV
325  bool colRed = false;
326 
327  m_low.clear();
328  h_crate_time_offsets->Set(0);
329 
330  for (unsigned short i = 0; i < 52; i++) {
331  std::string h_title = "ECL/time_crate_" + std::to_string(i + 1) + "_Thr1GeV";
332  h_time_crate_Thr1GeV = findHist(h_title);
333  if (h_time_crate_Thr1GeV != NULL) {
334  h_crate_time_offsets->SetPoint(i, i + 0.5, h_time_crate_Thr1GeV->GetMean());
335  h_crate_time_offsets->SetPointError(i, 0.1, h_time_crate_Thr1GeV->GetMeanError());
336  if (h_time_crate_Thr1GeV->Integral() > 100) {
337  double yval = (h_time_crate_Thr1GeV->GetMean() > 0) ?
338  h_time_crate_Thr1GeV->GetMean() - 2 * h_time_crate_Thr1GeV->GetMeanError() :
339  h_time_crate_Thr1GeV->GetMean() + 2 * h_time_crate_Thr1GeV->GetMeanError();
340  if (abs(yval) > m_CrateTimeOffsetsMax) colRed = true;
341  } else m_low.push_back(i + 1);
342  }
343  }
344 
345  c_crate_time_offsets->Clear();
346  c_crate_time_offsets->SetGrid();
347  c_crate_time_offsets->cd();
348 
349  c_crate_time_offsets->Pad()->SetFrameFillColor(10);
350  c_crate_time_offsets->Pad()->SetFillColor(kWhite);
351 
352  h_crate_time_offsets->SetMinimum(-50);
353  h_crate_time_offsets->SetMaximum(50);
354  h_crate_time_offsets->GetXaxis()->Set(52, 0, 52);
355  for (unsigned short i = 0; i < 52; i++) h_crate_time_offsets->GetXaxis()->SetBinLabel(i + 1, std::to_string(i + 1).c_str());
356  h_crate_time_offsets->GetXaxis()->LabelsOption("v");
357 
358  h_crate_time_offsets->Draw("AP");
359 
360  if (!m_low.empty()) {
361  auto tg = new TLatex(5, 40, "Low statistics");
362  tg->SetTextSize(.06);
363  tg->SetTextAlign(12);
364  tg->Draw();
365  c_crate_time_offsets->Pad()->SetFillColor(kYellow);
366  if (m_low.size() < 52) {
367  std::ostringstream sstream;
368  std::copy(m_low.begin(), m_low.end() - 1, std::ostream_iterator<short>(sstream, ","));
369  sstream << m_low.back();
370  std::string str = "Crates: " + sstream.str();
371  auto tg1 = new TLatex(5, 30, str.c_str());
372  tg1->SetTextSize(.03);
373  tg1->SetTextAlign(12);
374  tg1->Draw();
375  } else {
376  auto tg1 = new TLatex(5, 30, "All crates");
377  tg1->SetTextSize(.06);
378  tg1->SetTextAlign(12);
379  tg1->Draw();
380  }
381  }
382 
383  m_lower_boundary_time_offsets->Draw();
384  m_upper_boundary_time_offsets->Draw();
385 
386  if (colRed) c_crate_time_offsets->Pad()->SetFillColor(kRed);
387 
388  c_crate_time_offsets->Draw();
389  c_crate_time_offsets->Modified();
390  c_crate_time_offsets->Update();
391 
392  //DQM summary logic
393  h_logic_summary->Reset();
394 
395  TH1* h_fail_crateid = findHist("ECL/fail_crateid");
396  TH1* h_fail_shaperid = findHist("ECL/fail_shaperid");
397  if (h_fail_crateid != NULL && h_fail_shaperid != NULL) {
398  for (unsigned short i = 0; i < 52; i++) {
399  if (h_fail_crateid->GetBinContent(i + 1)) {
400  for (unsigned short j = i * 12; j < i * 12 + 12; j++) {
401  if (h_fail_shaperid->GetBinContent(j + 1)) {
402  unsigned short shaper = (j + 1) - 12 * i;
403  h_logic_summary->SetBinContent(h_logic_summary->FindBin(i + 1, shaper), h_fail_shaperid->GetBinContent(j + 1));
404  }
405  }
406  }
407  }
408  c_logic_summary->Clear();
409  c_logic_summary->SetGrid();
410  c_logic_summary->cd();
411 
412  c_logic_summary->Pad()->SetFrameFillColor(10);
413  c_logic_summary->Pad()->SetFillColor(kWhite);
414 
415  if (h_logic_summary->GetMaximum() > 0
416  && h_logic_summary->GetMaximum() < m_LogicTestMax) c_logic_summary->Pad()->SetFillColor(kYellow);
417  if (h_logic_summary->GetMaximum() >= m_LogicTestMax) c_logic_summary->Pad()->SetFillColor(kRed);
418  h_logic_summary->Draw("textcol");
419  c_logic_summary->Draw();
420  c_logic_summary->Modified();
421  c_logic_summary->Update();
422  }
423 }
424 
426 {
427  B2DEBUG(20, "DQMHistAnalysisECL: endRun called");
428 }
429 
430 
432 {
433  B2DEBUG(20, "terminate called");
434  delete c_quality_analysis;
436  delete c_bad_quality_analysis;
437  delete c_trigtag1_analysis;
438  delete c_trigtag2_analysis;
439  delete c_adc_hits_analysis;
445  delete m_lower_boundary_time_offsets;
446  delete m_upper_boundary_time_offsets;
447  delete c_crate_time_offsets;
448  delete h_time_crate_Thr1GeV;
449  delete h_crate_time_offsets;
450  delete c_logic_summary;
451  delete h_logic_summary;
452 }
Belle2::DQMHistAnalysisECLModule::beginRun
virtual void beginRun() override
Call when a run begins.
Definition: DQMHistAnalysisECL.cc:116
Belle2::DQMHistAnalysisECLModule::m_HitMapThresholds
std::vector< double > m_HitMapThresholds
Parameters for hit map histograms.
Definition: DQMHistAnalysisECL.h:74
Belle2::DQMHistAnalysisECLModule
This module is for analysis of ECL DQM histograms.
Definition: DQMHistAnalysisECL.h:50
Belle2::DQMHistAnalysisECLModule::c_quality_fit_data_analysis
TCanvas * c_quality_fit_data_analysis
TCanvas for quality_fit_data .
Definition: DQMHistAnalysisECL.h:125
Belle2::DQMHistAnalysisECLModule::m_lower_boundary_trigtag2
TLine * m_lower_boundary_trigtag2
TLine to show lower boundary for 'trigtag2_trigid' histogram.
Definition: DQMHistAnalysisECL.h:86
Belle2::DQMHistAnalysisECLModule::m_low
std::vector< short > m_low
Vector for crates IDs w/ low statistics.
Definition: DQMHistAnalysisECL.h:83
Belle2::DQMHistAnalysisECLModule::c_trigtag1_analysis
TCanvas * c_trigtag1_analysis
TCanvas for trigtag1 .
Definition: DQMHistAnalysisECL.h:115
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistAnalysisECLModule::m_WaveformOption
std::vector< std::string > m_WaveformOption
Options for waveform histograms.
Definition: DQMHistAnalysisECL.h:76
Belle2::DQMHistAnalysisECLModule::c_quality_analysis
TCanvas * c_quality_analysis
TCanvas for quality .
Definition: DQMHistAnalysisECL.h:109
Belle2::DQMHistAnalysisECLModule::c_quality_other_analysis
TCanvas * c_quality_other_analysis
TCanvas for quality_other .
Definition: DQMHistAnalysisECL.h:111
Belle2::DQMHistAnalysisECLModule::c_wf_analysis
std::vector< TCanvas * > c_wf_analysis
Vector of TCanvases for waveforms .
Definition: DQMHistAnalysisECL.h:131
Belle2::DQMHistAnalysisModule::findHist
static TH1 * findHist(const std::string &histname)
Find histogram.
Definition: DQMHistAnalysis.cc:83
Belle2::DQMHistAnalysisECLModule::h_time_crate_Thr1GeV
TH1 * h_time_crate_Thr1GeV
Histogram showing signal times from ECL crates (Thr.
Definition: DQMHistAnalysisECL.h:99
Belle2::DQMHistAnalysisECLModule::c_bad_quality_analysis
TCanvas * c_bad_quality_analysis
TCanvas for bad_quality .
Definition: DQMHistAnalysisECL.h:113
Belle2::DQMHistAnalysisECLModule::m_upper_boundary_trigtag2
TLine * m_upper_boundary_trigtag2
TLine to show upper boundary for 'trigtag2_trigid' histogram.
Definition: DQMHistAnalysisECL.h:88
Belle2::DQMHistAnalysisECLModule::terminate
virtual void terminate() override
Terminate.
Definition: DQMHistAnalysisECL.cc:431
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisECLModule::initialize
virtual void initialize() override
Initialize the module.
Definition: DQMHistAnalysisECL.cc:49
Belle2::DQMHistAnalysisECLModule::event
virtual void event() override
Event processor.
Definition: DQMHistAnalysisECL.cc:142
Belle2::DQMHistAnalysisECLModule::c_adc_hits_analysis
TCanvas * c_adc_hits_analysis
TCanvas for adc_hits .
Definition: DQMHistAnalysisECL.h:119
Belle2::DQMHistAnalysisECLModule::normalize
void normalize(TCanvas *, const std::string &, const Double_t &)
Normalize histograms.
Definition: DQMHistAnalysisECL.cc:121
Belle2::DQMHistAnalysisECLModule::c_cid_analysis
std::vector< TCanvas * > c_cid_analysis
Vector of TCanvases for hit map .
Definition: DQMHistAnalysisECL.h:128
Belle2::DQMHistAnalysisECLModule::c_crate_time_offsets
TCanvas * c_crate_time_offsets
TCanvas for time offsets.
Definition: DQMHistAnalysisECL.h:96
Belle2::DQMHistAnalysisECLModule::m_LogicTestMax
int m_LogicTestMax
Maximum of fails for logic test.
Definition: DQMHistAnalysisECL.h:80
Belle2::DQMHistAnalysisECLModule::m_CrateTimeOffsetsMax
double m_CrateTimeOffsetsMax
Maximum boundary for crate time offsets.
Definition: DQMHistAnalysisECL.h:78
Belle2::DQMHistAnalysisECLModule::c_ampfail_quality_analysis
TCanvas * c_ampfail_quality_analysis
TCanvas for ampfail_quality .
Definition: DQMHistAnalysisECL.h:121
Belle2::DQMHistAnalysisECLModule::c_logic_summary
TCanvas * c_logic_summary
TCanvas for ECL logic summary.
Definition: DQMHistAnalysisECL.h:104
Belle2::DQMHistAnalysisECLModule::h_crate_time_offsets
TGraphErrors * h_crate_time_offsets
Graph for time offsets.
Definition: DQMHistAnalysisECL.h:101
Belle2::DQMHistAnalysisECLModule::~DQMHistAnalysisECLModule
virtual ~DQMHistAnalysisECLModule()
Destructor.
Definition: DQMHistAnalysisECL.cc:47
Belle2::DQMHistAnalysisECLModule::h_logic_summary
TH2F * h_logic_summary
Histogram for ECL logic summary.
Definition: DQMHistAnalysisECL.h:106
Belle2::DQMHistAnalysisECLModule::c_timefail_quality_analysis
TCanvas * c_timefail_quality_analysis
TCanvas for timefail_quality .
Definition: DQMHistAnalysisECL.h:123
Belle2::DQMHistAnalysisECLModule::c_trigtag2_analysis
TCanvas * c_trigtag2_analysis
TCanvas for trigtag2 .
Definition: DQMHistAnalysisECL.h:117
Belle2::DQMHistAnalysisECLModule::endRun
virtual void endRun() override
Call when a run ends.
Definition: DQMHistAnalysisECL.cc:425
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27