Belle II Software  release-08-01-10
DQMHistAnalysisECLSummary.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/DQMHistAnalysisECLSummary.h>
11 
12 //ECL
13 #include <ecl/geometry/ECLGeometryPar.h>
14 
15 //ROOT
16 #include <TStyle.h>
17 #include <TLine.h>
18 #include <TText.h>
19 
20 //boost
21 #include "boost/format.hpp"
22 
23 using namespace Belle2;
24 
25 REG_MODULE(DQMHistAnalysisECLSummary);
26 
29  m_neighbours_obj("F", 0.1)
30 {
31 
32  B2DEBUG(20, "DQMHistAnalysisECLSummary: Constructor done.");
33  addParam("pvPrefix", m_pvPrefix, "Prefix to use for PVs registered by this module",
34  std::string("ECL:channels_info:"));
35  addParam("useChannelMask", m_useChannelMask,
36  "Mask Cell IDs based on information from ECL PVs",
37  true);
38  addParam("maxDeviationForOccupancy", m_maxDeviationForOccupancy,
39  "The higher this parameter, the larger differences in occupancy are allowed for adjacent channels", 0.28);
40  addParam("maxDeviationForChi2", m_maxDeviationForChi2,
41  "The higher this parameter, the larger differences in the number of hits with bad chi2 are allowed for adjacent channels",
42  2.5);
43  addParam("onlyIfUpdated", m_onlyIfUpdated, "If true (default), update EPICS PVs only if there histograms were updated.",
44  true);
45 }
46 
47 
49 {
50 }
51 
53 {
55 
56  //=== Set up ECL alarms and corresponding PVs
57 
58  m_ecl_alarms = {
59  {"dead", "#splitline{dead}{channels}", 1, 1, 1e5},
60  {"cold", "#splitline{cold}{channels}", 1, 2, 1e5},
61  {"hot", "#splitline{hot}{channels}", 25, 50, 1e5},
62  {"bad_chi2", "#splitline{bad #chi^{2}}{channels}", 5, 10, 1e6},
63  {"bad_fit", "#splitline{fit incon-}{sistencies}", 5, 10, 0 },
64  };
65 
66  // Prepare EPICS PVs
67  for (auto& alarm : m_ecl_alarms) {
68  // By crate
69  for (int crate_id = 1; crate_id <= ECL::ECL_CRATES; crate_id++) {
70  std::string pv_name = (boost::format("crate%02d:%s") % crate_id % alarm.name).str();
71  registerEpicsPV(m_pvPrefix + pv_name, pv_name, false);
72  }
73  // Totals
74  for (auto& ecl_part : {"All", "FWDEndcap", "Barrel", "BWDEndcap"}) {
75  std::string pv_name = (boost::format("%s:%s") % ecl_part % alarm.name).str();
76  registerEpicsPV(m_pvPrefix + pv_name, pv_name, false);
77  }
78  // Masked Cell IDs
79  std::string mask_pv_name = (boost::format("mask:%s") % alarm.name).str();
80  registerEpicsPV(m_pvPrefix + mask_pv_name, mask_pv_name, false);
81  }
82  updateEpicsPVs(5.0);
83 
85 
86  //=== Set up the histogram to indicate alarm status
87 
88  TString title = "#splitline{ECL errors monitoring}";
89  title += "{E - Error, W - Warning, L - Low statistics}";
90  title += ";ECLCollector ID (same as Crate ID)";
91  h_channels_summary = new TH2F("channels_summary", title,
92  ECL::ECL_CRATES, 1, ECL::ECL_CRATES + 1,
93  m_ecl_alarms.size(), 0, m_ecl_alarms.size());
94 
95  // Do not show statistics box.
96  h_channels_summary->SetStats(0);
97  h_channels_summary->SetMinimum(0);
98  h_channels_summary->SetMaximum(1);
99 
100  //=== Set X axis labels
101 
102  for (int i = 1; i <= ECL::ECL_CRATES; i++) {
103  h_channels_summary->GetXaxis()->SetBinLabel(i, std::to_string(i).c_str());
104  }
105  h_channels_summary->LabelsOption("v", "X"); // Rotate X axis labels 90 degrees
106  h_channels_summary->SetTickLength(0, "XY");
107 
108  //=== Customize offsets and margins
109 
110  h_channels_summary->GetXaxis()->SetTitleOffset(0.95);
111  h_channels_summary->GetXaxis()->SetTitleSize(0.05);
112  h_channels_summary->GetXaxis()->SetLabelSize(0.04);
113  h_channels_summary->GetYaxis()->SetLabelSize(0.06);
114 
115  c_channels_summary = new TCanvas("ECL/c_channels_summary_analysis");
116  c_channels_summary->SetTopMargin(0.10);
117  c_channels_summary->SetLeftMargin(0.20);
118  c_channels_summary->SetRightMargin(0.005);
119  c_channels_summary->SetBottomMargin(0.10);
120 
121  //=== Additional canvases/histograms to display which channels have problems
122  c_occupancy = new TCanvas("ECL/c_cid_Thr5MeV_overlaid_analysis");
123  c_bad_chi2 = new TCanvas("ECL/c_bad_quality_overlaid_analysis");
124 
125  h_bad_occ_overlay = new TH1F("bad_occ_overlay", "", ECL::ECL_TOTAL_CHANNELS,
126  1, ECL::ECL_TOTAL_CHANNELS + 1);
127  h_bad_occ_overlay->SetLineColor(kRed);
128  h_bad_occ_overlay->SetLineStyle(kDashed);
129  h_bad_occ_overlay->SetFillColor(kRed);
130  h_bad_occ_overlay->SetFillStyle(3007);
131 
132  h_bad_chi2_overlay = new TH1F("bad_chi2_overlay", "", ECL::ECL_TOTAL_CHANNELS,
133  1, ECL::ECL_TOTAL_CHANNELS + 1);
134  h_bad_chi2_overlay->SetLineColor(kRed);
135  h_bad_chi2_overlay->SetLineStyle(kDashed);
136  h_bad_chi2_overlay->SetFillColor(kRed);
137  h_bad_chi2_overlay->SetFillStyle(3007);
138 
139  h_bad_occ_overlay_green = new TH1F("bad_occ_overlay_green", "", ECL::ECL_TOTAL_CHANNELS,
140  1, ECL::ECL_TOTAL_CHANNELS + 1);
141  h_bad_occ_overlay_green->SetLineColor(kGreen);
142  h_bad_occ_overlay_green->SetLineStyle(kDotted);
143  h_bad_occ_overlay_green->SetFillColor(kGreen);
144  h_bad_occ_overlay_green->SetFillStyle(3013);
145 
146  h_bad_chi2_overlay_green = new TH1F("bad_chi2_overlay_green", "", ECL::ECL_TOTAL_CHANNELS,
147  1, ECL::ECL_TOTAL_CHANNELS + 1);
148  h_bad_chi2_overlay_green->SetLineColor(kGreen);
149  h_bad_chi2_overlay_green->SetLineStyle(kDotted);
150  h_bad_chi2_overlay_green->SetFillColor(kGreen);
151  h_bad_chi2_overlay_green->SetFillStyle(3013);
152 
153  B2DEBUG(20, "DQMHistAnalysisECLSummary: initialized.");
154 }
155 
157 {
158  B2DEBUG(20, "DQMHistAnalysisECLSummary: beginRun called.");
159 
160  //=== Update m_ecl_alarms based on PV limits
161 
163 
164  //=== Set Y axis labels from m_ecl_alarms
165 
166  for (size_t i = 0; i < m_ecl_alarms.size(); i++) {
167  TString label = m_ecl_alarms[i].title;
168  label += " #geq ";
169  label += m_ecl_alarms[i].alarm_limit;
170  h_channels_summary->GetYaxis()->SetBinLabel(i + 1, label);
171  }
172 }
173 
175 {
176  h_channels_summary->Reset();
177 
178  TH1* h_total_events = findHist("ECL/event", m_onlyIfUpdated);
179  if (!h_total_events) return;
180  m_total_events = h_total_events->GetEntries();
181 
182  // [alarm_type][crate_id - 1] -> number of problematic channels in that crate
183  std::vector< std::vector<int> > alarm_counts = updateAlarmCounts();
184 
185  //=== Set warning and error indicators on the histogram
186 
187  const double HISTCOLOR_RED = 0.9;
188  const double HISTCOLOR_GREEN = 0.45;
189  const double HISTCOLOR_ORANGE = 0.65;
190  const double HISTCOLOR_BLUE = 0.01;
191 
192  std::vector<TText*> labels;
193 
194  for (size_t alarm_idx = 0; alarm_idx < alarm_counts.size(); alarm_idx++) {
195  for (size_t crate = 0; crate < alarm_counts[alarm_idx].size(); crate++) {
196  double color;
197  char label_text[2] = {0};
198  int alarm_limit = m_ecl_alarms[alarm_idx].alarm_limit;
199  int warning_limit = m_ecl_alarms[alarm_idx].warning_limit;
200  if (m_total_events < m_ecl_alarms[alarm_idx].required_statistics) {
201  color = HISTCOLOR_BLUE;
202  label_text[0] = 'L';
203  } else if (alarm_counts[alarm_idx][crate] >= alarm_limit) {
204  color = HISTCOLOR_RED;
205  label_text[0] = 'E';
206  } else if (alarm_counts[alarm_idx][crate] >= warning_limit) {
207  color = HISTCOLOR_ORANGE;
208  label_text[0] = 'W';
209  } else {
210  color = HISTCOLOR_GREEN;
211  }
212  if (label_text[0] == 'E' || label_text[0] == 'W') {
213  B2DEBUG(100, "Non-zero (" << alarm_counts[alarm_idx][crate]
214  << ") for alarm_idx, crate = " << alarm_idx << ", " << crate);
215  }
216  h_channels_summary->SetBinContent(crate + 1, alarm_idx + 1, color);
217  if (label_text[0] != 0) {
218  auto text = new TText((crate + 1.5), (alarm_idx + 0.5), label_text);
219  text->SetTextColor(kWhite);
220  text->SetTextSize(0.03);
221  text->SetTextAlign(22); // centered
222  labels.push_back(text);
223  }
224  }
225  }
226 
227  //=== Draw histogram, labels and grid
228 
229  // Customize title
230  int gstyle_title_h = gStyle->GetTitleH();
231  int gstyle_title_x = gStyle->GetTitleX();
232  int gstyle_title_y = gStyle->GetTitleY();
233  gStyle->SetTitleH(0.04);
234  gStyle->SetTitleX(0.60);
235  gStyle->SetTitleY(1.00);
236 
237  c_channels_summary->Clear();
238  c_channels_summary->cd();
239 
240  //=== Prepare special style objects to use correct color palette
241  // and use it only for this histogram
242  // Based on https://root-forum.cern.ch/t/different-color-palettes-for-different-plots-with-texec/5250/3
243  // and https://root.cern/doc/master/multipalette_8C.html
244 
245  if (!m_ecl_style) delete m_ecl_style;
246  if (!m_default_style) delete m_default_style;
247 
248  m_ecl_style = new TExec("ecl_style",
249  "gStyle->SetPalette(kRainBow);"
250  "channels_summary->SetDrawOption(\"col\");");
251  h_channels_summary->GetListOfFunctions()->Add(m_ecl_style);
252  m_default_style = new TExec("default_style",
253  "gStyle->SetPalette(kBird);");
254 
255  //=== Draw with special style
256  // https://root.cern.ch/js/latest/examples.htm#th2_colpal77
257  h_channels_summary->Draw("");
258  h_channels_summary->Draw("colpal55;same");
259  for (auto& text : labels) {
260  text->Draw();
261  }
263  m_default_style->Draw("same");
264 
265  //
266  c_channels_summary->Modified();
267  c_channels_summary->Update();
268  c_channels_summary->Draw();
269 
270  gStyle->SetTitleH(gstyle_title_h);
271  gStyle->SetTitleX(gstyle_title_x);
272  gStyle->SetTitleY(gstyle_title_y);
273 }
274 
276 {
277  B2DEBUG(20, "DQMHistAnalysisECLSummary: endRun called");
278 
279  updateAlarmCounts(true);
280 }
281 
282 
284 {
285  B2DEBUG(20, "terminate called");
286  delete c_channels_summary;
287  delete c_occupancy;
288  delete c_bad_chi2;
289  delete h_bad_occ_overlay;
290  delete h_bad_chi2_overlay;
291 }
292 
293 std::pair<int, DQMHistAnalysisECLSummaryModule::ECLAlarmType> DQMHistAnalysisECLSummaryModule::getAlarmByName(std::string name)
294 {
295  int index = 0;
296  for (auto& alarm_info : m_ecl_alarms) {
297  if (alarm_info.name == name) return {index, alarm_info};
298  index++;
299  }
300  B2FATAL("Could not get ECL alarm " + name);
301  return {-1, m_ecl_alarms[0]};
302 }
303 
305 {
306  //===== Get alarm limits
307 
308  for (auto& alarm : m_ecl_alarms) {
309  // In the current version, use only first crate PV to get alarm limits
310  int crate_id = 1;
311 
312  std::string pv_name = (boost::format("crate%02d:%s") % crate_id % alarm.name).str();
313  double unused = 0;
314  double upper_warning_limit = 0;
315  double upper_alarm_limit = 0;
316  bool accessed = requestLimitsFromEpicsPVs(pv_name, unused, unused, upper_warning_limit, upper_alarm_limit);
317 
318  if (!accessed || upper_alarm_limit <= 0 || upper_warning_limit <= 0) {
319  B2WARNING("Failed to get alarm limits");
320  return;
321  }
322 
323  alarm.alarm_limit = upper_alarm_limit;
324  alarm.warning_limit = upper_warning_limit;
325  }
326 
327  //===== Get masked channels
328 
329  // alarm_type -> array_of_masked_channels
330  static std::map<std::string, dbr_sts_long_array> mask_info;
331  if (!getMaskedChannels(mask_info)) {
332  B2WARNING("Failed to get arrays of masked channels");
333  return;
334  }
335 
336  for (auto& alarm : m_ecl_alarms) {
337  // Update the list of masked channels
338  m_mask[alarm.name].clear();
339  for (int i = 0; i < c_max_masked_channels; i++) {
340  int cell_id = mask_info[alarm.name].value[i];
341  // The array of masked channels goes like this: [15, 22, 716, 0, 0, 0, 0, 0]
342  // So if 0 (invalid channel ID) is encountered, stop reading further.
343  if (cell_id == 0) break;
344  m_mask[alarm.name].insert(cell_id);
345  }
346  }
347 }
348 
349 bool DQMHistAnalysisECLSummaryModule::getMaskedChannels(std::map<std::string, dbr_sts_long_array>& mask_info)
350 {
351  for (auto& alarm : m_ecl_alarms) {
352  std::string mask_pv_name = (boost::format("mask:%s") % alarm.name).str();
353  chid mask_chid = getEpicsPVChID(mask_pv_name);
354 
355  if (mask_chid == nullptr) return false;
356  // ca_array_get data is NOT VALID here.
357  // It will only be valid after ca_pend_io returns ECA_NORMAL
358  auto r = ca_array_get(DBR_STS_LONG, c_max_masked_channels, mask_chid, &mask_info[alarm.name]);
359  if (r != ECA_NORMAL) return false;
360  }
361 
362  // ca_pend_io has to be called here,
363  // see https://epics.anl.gov/base/R3-14/12-docs/CAref.html#ca_get
364 
365  auto r = ca_pend_io(5.0);
366  if (r != ECA_NORMAL) return false;
367 
368  return true;
369 }
370 
371 std::vector< std::vector<int> > DQMHistAnalysisECLSummaryModule::updateAlarmCounts(bool update_mirabelle)
372 {
373  std::vector< std::vector<int> > alarm_counts(
374  m_ecl_alarms.size(), std::vector<int>(ECL::ECL_CRATES));
375 
376  //=== Get number of fit inconsistencies
377 
378  TH1* h_fail_crateid = findHist("ECL/fail_crateid", m_onlyIfUpdated);
379 
380  const int fit_alarm_index = getAlarmByName("bad_fit").first;
381  for (int crate_id = 1; crate_id <= ECL::ECL_CRATES; crate_id++) {
382  int errors_count = 0;
383  if (h_fail_crateid) {
384  errors_count = h_fail_crateid->GetBinContent(crate_id);
385  } else {
386  errors_count = 9999;
387  }
388 
389  alarm_counts[fit_alarm_index][crate_id - 1] += errors_count;
390  }
391 
392  // [Cell ID] -> error_bitmask
393  std::map<int, int> error_bitmasks;
394 
395  //=== Get number of dead/cold/hot channels
396  // cppcheck-suppress unassignedVariable
397  for (auto& [cell_id, error_bitmask] : getChannelsWithOccupancyProblems()) {
398  error_bitmasks[cell_id] |= error_bitmask;
399  }
400  //=== Get number of channels with bad_chi2
401  // cppcheck-suppress unassignedVariable
402  for (auto& [cell_id, error_bitmask] : getChannelsWithChi2Problems()) {
403  error_bitmasks[cell_id] |= error_bitmask;
404  }
405 
406  //=== Combine the information
407 
408  if (!update_mirabelle) {
409  h_bad_chi2_overlay->Reset();
410  h_bad_occ_overlay->Reset();
411  h_bad_chi2_overlay_green->Reset();
412  h_bad_occ_overlay_green->Reset();
413  }
414 
415  static std::vector<std::string> indices = {"dead", "cold", "hot", "bad_chi2"};
416  for (auto& index_name : indices) {
417  int alarm_index = getAlarmByName(index_name).first;
418  int alarm_bit = 1 << alarm_index;
419  for (auto& [cid, error_bitmask] : error_bitmasks) {
420  int crate_id = m_mapper.getCrateID(cid);
421  if ((error_bitmask & alarm_bit) == 0) continue;
422 
423  bool masked = (m_mask[index_name].find(cid) != m_mask[index_name].end());
424 
425  if (!masked) alarm_counts[alarm_index][crate_id - 1] += 1;
426 
427  if (update_mirabelle) continue;
428 
429  if (index_name == "bad_chi2") {
430  if (!masked) h_bad_chi2_overlay->SetBinContent(cid, 1);
431  else h_bad_chi2_overlay_green->SetBinContent(cid, 1);
432  } else if (index_name == "dead" || index_name == "cold" || index_name == "hot") {
433  if (!masked) h_bad_occ_overlay->SetBinContent(cid, 1);
434  else h_bad_occ_overlay_green->SetBinContent(cid, 1);
435  }
436  }
437 
438  if (update_mirabelle) continue;
439 
440  //=== Prepare the histogram to display the list of all bad channels
441  if (index_name == "hot" || index_name == "bad_chi2") {
442  TCanvas* current_canvas;
443  TH1* main_hist;
444  TH1F* overlay_hist;
445  TH1F* overlay_hist_green;
446  if (index_name == "hot") {
447  main_hist = findHist("ECL/cid_Thr5MeV", m_onlyIfUpdated);
448  overlay_hist = h_bad_occ_overlay;
449  overlay_hist_green = h_bad_occ_overlay_green;
450  current_canvas = c_occupancy;
451  } else {
452  main_hist = findHist("ECL/bad_quality", m_onlyIfUpdated);
453  overlay_hist = h_bad_chi2_overlay;
454  overlay_hist_green = h_bad_chi2_overlay_green;
455  current_canvas = c_bad_chi2;
456  }
457 
458  if (main_hist) {
459  for (auto& overlay : {overlay_hist, overlay_hist_green}) {
460  for (int bin_id = 1; bin_id <= ECL::ECL_TOTAL_CHANNELS; bin_id++) {
461  if (overlay->GetBinContent(bin_id) == 0) continue;
462  // Do not adjust bin height for dead channels
463  if (main_hist->GetBinContent(bin_id) == 0) continue;
464  overlay->SetBinContent(bin_id, main_hist->GetBinContent(bin_id));
465  }
466  }
467 
468  current_canvas->Clear();
469  current_canvas->cd();
470  main_hist->Draw("hist");
471  overlay_hist->Draw("hist;same");
472  overlay_hist_green->Draw("hist;same");
473  current_canvas->Modified();
474  current_canvas->Update();
475  current_canvas->Draw();
476  }
477  }
478  }
479 
480  //== Update EPICS PVs or MiraBelle monObjs
481 
482  for (size_t alarm_idx = 0; alarm_idx < alarm_counts.size(); alarm_idx++) {
483  auto& alarm = m_ecl_alarms[alarm_idx];
484  std::map<std::string, int> total;
485  // Convert values per crate to totals
486  for (size_t crate = 0; crate < alarm_counts[alarm_idx].size(); crate++) {
487  int crate_id = crate + 1;
488  int value = alarm_counts[alarm_idx][crate];
489 
490  if (!update_mirabelle) {
491  std::string pv_name = (boost::format("crate%02d:%s") % crate_id % alarm.name).str();
492  setEpicsPV(pv_name, value);
493  }
494 
495  total["All"] += value;
496  if (crate_id <= ECL::ECL_BARREL_CRATES) {
497  total["Barrel"] += value;
498  } else if (crate_id <= ECL::ECL_BARREL_CRATES + ECL::ECL_FWD_CRATES) {
499  total["FWDEndcap"] += value;
500  } else {
501  total["BWDEndcap"] += value;
502  }
503  }
504  // Export totals
505  for (auto& ecl_part : {"All", "FWDEndcap", "Barrel", "BWDEndcap"}) {
506  std::string pv_name = (boost::format("%s:%s") % ecl_part % alarm.name).str();
507  if (update_mirabelle) {
508  std::string var_name = pv_name;
509  std::replace(var_name.begin(), var_name.end(), ':', '_');
510  B2DEBUG(100, var_name << " = " << total[ecl_part]);
511  m_monObj->setVariable(var_name, total[ecl_part]);
512  } else {
513  setEpicsPV(pv_name, total[ecl_part]);
514  }
515  }
516  }
517 
518  if (!update_mirabelle) updateEpicsPVs(5.0);
519 
520  return alarm_counts;
521 }
522 
524 {
525  static std::vector< std::vector<short> > neighbours(ECL::ECL_TOTAL_CHANNELS);
526  if (neighbours[0].size() == 0) {
527  for (int cid0 = 0; cid0 < ECL::ECL_TOTAL_CHANNELS; cid0++) {
528  // [0]First is the crystal itself.
529  // [1]Second is PHI neighbour, phi+1.
530  // [2]Third is PHI neighbour, phi-1.
531  // [3]Next one (sometimes two) are THETA neighbours in theta-1.
532  // [4]Next one (sometimes two) are THETA neighbours in theta+1.
533  neighbours[cid0] = m_neighbours_obj.getNeighbours(cid0 + 1);
534  // Remove first element (the crystal itself)
535  neighbours[cid0].erase(neighbours[cid0].begin());
536  }
537  }
538 
539  TH1* h_occupancy = findHist("ECL/cid_Thr5MeV", m_onlyIfUpdated);
540  const double max_deviation = m_maxDeviationForOccupancy;
541  return getSuspiciousChannels(h_occupancy, m_total_events, neighbours,
542  max_deviation, true);
543 }
544 
546 {
548 
549  static std::vector< std::vector<short> > neighbours(ECL::ECL_TOTAL_CHANNELS);
550  if (neighbours[0].size() == 0) {
551  for (int cid_center = 1; cid_center <= ECL::ECL_TOTAL_CHANNELS; cid_center++) {
552  geom->Mapping(cid_center - 1);
553  const int theta_id_center = geom->GetThetaID();
554  int phi_id_center = geom->GetPhiID();
555  const int crystals_in_ring = m_neighbours_obj.getCrystalsPerRing(theta_id_center);
556  phi_id_center = phi_id_center * 144 / crystals_in_ring;
557  for (int cid0 = 0; cid0 < 8736; cid0++) {
558  if (cid0 == cid_center - 1) continue;
559  geom->Mapping(cid0);
560  int theta_id = geom->GetThetaID();
561  int phi_id = geom->GetPhiID();
562  phi_id = phi_id * 144 / m_neighbours_obj.getCrystalsPerRing(theta_id);
563  if (std::abs(theta_id - theta_id_center) <= 2 &&
564  std::abs(phi_id - phi_id_center) <= 2) {
565  neighbours[cid_center - 1].push_back(cid0);
566  }
567  }
568  }
569  }
570 
571  TH1* h_bad_chi2 = findHist("ECL/bad_quality", m_onlyIfUpdated);
572  const double max_deviation = m_maxDeviationForChi2;
573  return getSuspiciousChannels(h_bad_chi2, m_total_events, neighbours,
574  max_deviation, false);
575 }
576 
578  TH1* hist, double total_events,
579  const std::vector< std::vector<short> >& neighbours,
580  double max_deviation, bool occupancy_histogram)
581 {
582  std::map<int, int> retval;
583 
584  if (!hist) return retval;
585  //== A bit hacky solution to skip incorrectly
586  // filled histograms
587  if (hist->Integral() <= 0) return retval;
588 
589  //=== Extract alarm details
590  // cppcheck-suppress unassignedVariable
591  const auto& [dead_index, dead_alarm] = getAlarmByName("dead");
592  // cppcheck-suppress unassignedVariable
593  const auto& [cold_index, cold_alarm] = getAlarmByName("cold");
594  // cppcheck-suppress unassignedVariable
595  const auto& [hot_index, hot_alarm ] = getAlarmByName("hot");
596  // cppcheck-suppress unassignedVariable
597  const auto& [chi2_index, chi2_alarm] = getAlarmByName("bad_chi2");
598 
599  double min_required_events;
600 
601  if (occupancy_histogram) {
602  min_required_events = std::min({
603  dead_alarm.required_statistics,
604  cold_alarm.required_statistics,
605  hot_alarm.required_statistics,
606  });
607  } else {
608  min_required_events = chi2_alarm.required_statistics;
609  }
610 
611  if (total_events < min_required_events) return retval;
612 
613  int dead_bit = 1 << dead_index;
614  int cold_bit = 1 << cold_index;
615  int hot_bit = 1 << hot_index;
616  int chi2_bit = 1 << chi2_index;
617 
618  // == Search for dead channels
619 
620  if (occupancy_histogram) {
621  // This indicates that DQMHistAnalysisECL module is not included in the path
622  bool not_normalized = (findCanvas("ECL/c_cid_Thr5MeV_analysis") == nullptr);
623  if (total_events >= dead_alarm.required_statistics) {
624  // For null runs, it should be higher than 0.0001%
625  // For cosmic runs, number of hits with E > threshold should be higher than 0.01%
626  // (for physics runs, as opposed to cosmics, this can actually be set to higher value)
627  double min_occupancy = 1e-6;
628  if (not_normalized) {
629  // The histogram is not normalized, multiply the threshold by evt count
630  min_occupancy *= total_events;
631  }
632  for (int cid = 1; cid <= ECL::ECL_TOTAL_CHANNELS; cid++) {
633  if (hist->GetBinContent(cid) > min_occupancy) continue;
634  retval[cid] |= dead_bit;
635  }
636  }
637  if (total_events >= hot_alarm.required_statistics) {
638  // Number of hits with E > threshold should be less than 70%
639  double max_occupancy = 0.7;
640  if (not_normalized) {
641  // The histogram is not normalized, multiply the threshold by evt count
642  max_occupancy *= total_events;
643  }
644  for (int cid = 1; cid <= ECL::ECL_TOTAL_CHANNELS; cid++) {
645  if (hist->GetBinContent(cid) < max_occupancy) continue;
646  retval[cid] |= hot_bit;
647  }
648  }
649  }
650 
651  // == Search for cold and hot channels (or channels with bad chi2)
652 
653  for (int cid = 1; cid <= ECL::ECL_TOTAL_CHANNELS; cid++) {
654  double actual_value = hist->GetBinContent(cid);
655 
656  std::vector<short> neighb = neighbours[cid - 1];
657  std::multiset<double> values_sorted;
658  for (auto& neighbour_cid : neighb) {
659  values_sorted.insert(hist->GetBinContent(neighbour_cid));
660  }
661 
662  double expected_value;
663  if (occupancy_histogram) {
664  // Use "lower 25%" value:
665  expected_value = *std::next(values_sorted.begin(), 1 * neighb.size() / 4);
666  } else {
667  // Use "upper ~75%" value (48 is the expected number of neighbours)
668  expected_value = *std::next(values_sorted.begin(), 35 * neighb.size() / 48);
669  }
670  double deviation = std::abs((actual_value - expected_value) /
671  expected_value);
672 
673  if (!occupancy_histogram) {
674  // Min occupancy for high-energy (> 1GeV) hits
675  double min_occupancy = 1.41e-5;
676  if (findCanvas("ECL/c_bad_quality_analysis") == nullptr) {
677  // The histogram is not normalized, multiply the threshold by evt count
678  min_occupancy *= total_events;
679  }
680 
681  if (actual_value < min_occupancy) continue;
682  }
683 
684  if (deviation < max_deviation) continue;
685 
686  if (occupancy_histogram) {
687  if (actual_value < expected_value) retval[cid] |= cold_bit;
688  if (actual_value > 2 * expected_value) retval[cid] |= hot_bit;
689  } else {
690  if (actual_value > expected_value) retval[cid] |= chi2_bit;
691  }
692  }
693 
694  return retval;
695 }
696 
698 {
699  int x_min = hist->GetXaxis()->GetXmin();
700  int x_max = hist->GetXaxis()->GetXmax();
701  int y_min = hist->GetYaxis()->GetXmin();
702  int y_max = hist->GetYaxis()->GetXmax();
703  for (int x = x_min + 1; x < x_max; x++) {
704  auto l = new TLine(x, 0, x, 5);
705  l->SetLineStyle(kDashed);
706  l->Draw();
707  }
708  for (int y = y_min + 1; y < y_max; y++) {
709  auto l = new TLine(1, y, ECL::ECL_CRATES + 1, y);
710  l->SetLineStyle(kDashed);
711  l->Draw();
712  }
713 }
714 
std::map< int, int > getSuspiciousChannels(TH1 *hist, double total_events, const std::vector< std::vector< short > > &neighbours, double max_deviation, bool occupancy_histogram)
Get outlier channels, ones with values (occupancy, errors, etc) that are much higher than the values ...
double m_total_events
Number of events with ECL data in the current run.
double m_maxDeviationForChi2
The higher this parameter, the larger differences in the number of hits with bad chi2 are allowed for...
TH1F * h_bad_occ_overlay
Overlay to indicate bad channels on occupancy histogram.
void initialize() override final
Initialize the module.
void updateAlarmConfig()
Set alarm limits in DQM based on EPICS PV limits.
std::vector< std::vector< int > > updateAlarmCounts(bool update_mirabelle=false)
Calculate and return number of problems per crate for all alarm categories.
ECL::ECLChannelMapper m_mapper
Object to map ECL crystal ID to ECL crate ID.
void drawGrid(TH2 *hist)
Draw grid with TLine primitives for the specified histogram.
TCanvas * c_channels_summary
TCanvas for ECL alarms regarding suspicious channels.
static const int c_max_masked_channels
Size of an array with masked channels.
std::pair< int, ECLAlarmType > getAlarmByName(std::string name)
Returns index of the specific alarm type and detailed information.
TH2F * h_channels_summary
Summarized information about ECL channels.
double m_maxDeviationForOccupancy
The higher this parameter, the larger differences in occupancy are allowed for adjacent channels.
ECL::ECLNeighbours m_neighbours_obj
Object to get ECL crystal neighbours.
TCanvas * c_occupancy
ECL occupancy histogram with highlighted suspicious channels.
TExec * m_ecl_style
Special object to specify style parameters.
std::string m_pvPrefix
Prefix to use for PVs registered by this module.
MonitoringObject * m_monObj
MiraBelle monitoring object.
void event() override final
Event processor.
TCanvas * c_bad_chi2
ECL bad chi2 histogram with highlighted suspicious channels.
TH1F * h_bad_chi2_overlay
Overlay to indicate bad channels on chi2 histogram.
void endRun() override final
Call when a run ends.
bool m_onlyIfUpdated
If true (default), update EPICS PVs only if there were changes in the histograms.
TH1F * h_bad_chi2_overlay_green
Overlay to indicate masked bad channels on chi2 histogram.
void beginRun() override final
Call when a run begins.
DQMHistAnalysisECLSummaryModule()
< derived from DQMHistAnalysisModule class.
TExec * m_default_style
Special object to revert changes done by ecl_style.
TH1F * h_bad_occ_overlay_green
Overlay to indicate masked bad channels on occupancy histogram.
std::map< std::string, std::set< int > > m_mask
List of masked Cell IDs for each alarm type.
bool getMaskedChannels(std::map< std::string, dbr_sts_long_array > &mask_info)
Get the array of masked channels for each type of alarm case monitored by this module.
bool m_useChannelMask
If true, mask Cell IDs based on DQM:ECL:channels_info:{alarm_type} PV info.
std::vector< ECLAlarmType > m_ecl_alarms
Different alarms monitored in h_channels_summary.
The base class for the histogram analysis module.
TCanvas * findCanvas(TString cname)
Find canvas by name.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
chid getEpicsPVChID(std::string keyname)
Get EPICS PV Channel Id.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
int updateEpicsPVs(float timeout)
Update all EPICS PV (flush to network)
bool initFromFile()
Initialize channel mapper using data stored in default location.
int getCrateID(int iCOPPERNode, int iFINESSE, bool pcie40=false)
Get crate number by given COPPER node number and FINESSE number.
The Class for ECL Geometry Parameters.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
const std::vector< short int > & getNeighbours(short int cid) const
Return the neighbours for a given cell ID.
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Definition: ECLNeighbours.h:39
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.