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