Belle II Software  release-08-01-10
DQMHistAnalysisKLM2.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 /* Own header. */
10 #include <dqm/analysis/modules/DQMHistAnalysisKLM2.h>
11 
12 /* Basf2 headers. */
13 #include <klm/dataobjects/KLMChannelIndex.h>
14 
15 /* ROOT headers. */
16 #include <TROOT.h>
17 #include <TStyle.h>
18 #include <TColor.h>
19 
20 /* Standard Library*/
21 #include <algorithm>
22 
23 using namespace Belle2;
24 
25 REG_MODULE(DQMHistAnalysisKLM2);
26 
29  m_EklmElementNumbers{&(EKLMElementNumbers::Instance())}
30 {
31  setDescription("Module used to analyze KLM Efficiency DQM histograms (depends on tracking variables).");
32  addParam("HistogramDirectoryName", m_histogramDirectoryName, "Name of histogram directory", std::string("KLMEfficiencyDQM"));
33  addParam("MinEvents", m_minEvents, "Minimum events for delta histogram update", 5000000.);
34  addParam("RefHistoFile", m_refFileName, "Reference histogram file name", std::string("KLM_DQM_REF_BEAM.root"));
35  addParam("AlarmThreshold", m_alarmThr, "Set alarm threshold", float(0.9));
36  addParam("WarnThreshold", m_warnThr, "Set warn threshold", float(0.92));
37  addParam("Min2DEff", m_min, "2D efficiency min", float(0.5));
38  addParam("Max2DEff", m_max, "2D efficiency max", float(2));
39  addParam("RatioPlot", m_ratio, "2D efficiency ratio or difference plot ", bool(true));
40 
41  m_PlaneLine.SetLineColor(kMagenta);
42  m_PlaneLine.SetLineWidth(1);
43  m_PlaneLine.SetLineStyle(2); // dashed
44  m_PlaneText.SetTextAlign(22); // centered, middle
45  m_PlaneText.SetTextColor(kMagenta);
46  m_PlaneText.SetTextFont(42); // Helvetica regular
47  m_PlaneText.SetTextSize(0.02); // 2% of TPad's full height
48 }
49 
50 
51 
53 {
55 
56  //register EPICS PVs
57  registerEpicsPV("KLMEff:nEffBKLMLayers", "nEffBKLMLayers");
58  registerEpicsPV("KLMEff:nEffEKLMLayers", "nEffEKLMLayers");
59  registerEpicsPV("KLMEff:2DEffSettings", "2DEffSettings");
60  updateEpicsPVs(5.0);
61 
62  if (m_refFileName != "") {
63  m_refFile = TFile::Open(m_refFileName.data(), "READ");
64  }
65 
66  // Get Reference Histograms
67  if (m_refFile && m_refFile->IsOpen()) {
68  B2INFO("DQMHistAnalysisKLM2: reference root file (" << m_refFileName << ") FOUND, able to read ref histograms");
69 
70  m_ref_efficiencies_bklm = (TH1F*)m_refFile->Get((m_histogramDirectoryName + "/eff_bklm_plane").data());
71  m_ref_efficiencies_bklm->SetLineColor(2);
72  m_ref_efficiencies_bklm->SetOption("HIST");
73  m_ref_efficiencies_bklm->SetStats(false);
74 
75  m_ref_efficiencies_eklm = (TH1F*)m_refFile->Get((m_histogramDirectoryName + "/eff_eklm_plane").data());
76  m_ref_efficiencies_eklm->SetLineColor(2);
77  m_ref_efficiencies_eklm->SetOption("HIST");
78  m_ref_efficiencies_eklm->SetStats(false);
79 
80  } else {
81  B2WARNING("DQMHistAnalysisKLM2: reference root file (" << m_refFileName << ") not found, or closed");
82 
83  // Switch to absolute 2D efficiencies if reference histogram is not found
84  m_alarmThr = 0;
85  m_warnThr = 0;
86  m_ref_efficiencies_bklm = new TH1F("eff_bklm_plane", "Plane Efficiency in BKLM", BKLMElementNumbers::getMaximalLayerGlobalNumber(),
88  for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerGlobalNumber(); lay_id++) {
89  if (m_ratio) {
90  m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 1);
91  } else {
92  m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 0);
93  }
94  }
95 
96  m_ref_efficiencies_eklm = new TH1F("eff_eklm_plane", "Plane Efficiency in EKLM", EKLMElementNumbers::getMaximalPlaneGlobalNumber(),
98  for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalPlaneGlobalNumber(); lay_id++) {
99  if (m_ratio) {
100  m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 1);
101  } else {
102  m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 0);
103  }
104  }
105 
106  }
107  gROOT->cd();
108  m_c_eff_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_plane").data());
109  m_c_eff_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_plane").data());
110  m_c_eff_bklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_sector").data());
111  m_c_eff_eklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_sector").data());
112 
113  m_c_eff2d_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_bklm").data());
114  m_c_eff2d_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_eklm").data());
115 
116  m_eff_bklm = new TH1F((m_histogramDirectoryName + "/eff_bklm_plane").data(),
117  "Plane Efficiency in BKLM",
119  m_eff_bklm->GetXaxis()->SetTitle("Layer number");
120  m_eff_bklm->SetStats(false);
121  m_eff_bklm->SetOption("HIST");
122 
123  m_eff_eklm = new TH1F((m_histogramDirectoryName + "/eff_eklm_plane").data(),
124  "Plane Efficiency in EKLM",
126  m_eff_eklm->GetXaxis()->SetTitle("Plane number");
127  m_eff_eklm->SetStats(false);
128  m_eff_eklm->SetOption("HIST");
129 
130  m_eff_bklm_sector = new TH1F((m_histogramDirectoryName + "/eff_bklm_sector").data(),
131  "Sector Efficiency in BKLM",
133  m_eff_bklm_sector->GetXaxis()->SetTitle("Sector number");
134  m_eff_bklm_sector->SetStats(false);
135  m_eff_bklm_sector->SetOption("HIST");
136 
137  m_eff_eklm_sector = new TH1F((m_histogramDirectoryName + "/eff_eklm_sector").data(),
138  "Sector Efficiency in EKLM",
140  m_eff_eklm_sector->GetXaxis()->SetTitle("Sector number");
141  m_eff_eklm_sector->SetStats(false);
142  m_eff_eklm_sector->SetOption("HIST");
143 
144  /* register plots for delta histogramming */
145  // all ext hits
146  addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLM", HistDelta::c_Events, m_minEvents, 1);
147  addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLM", HistDelta::c_Events, m_minEvents, 1);
148  addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLMSector", HistDelta::c_Events, m_minEvents, 1);
149  addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLMSector", HistDelta::c_Events, m_minEvents, 1);
150 
151  // matched hits
152  addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLM", HistDelta::c_Events, m_minEvents, 1);
153  addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLM", HistDelta::c_Events, m_minEvents, 1);
154  addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLMSector", HistDelta::c_Events, m_minEvents, 1);
155  addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLMSector", HistDelta::c_Events, m_minEvents, 1);
156 
157  // 2D Efficiency Histograms
158  TString eff2d_hist_bklm_title;
159  TString eff2d_hist_eklm_title;
160  if (m_ratio) {
161  eff2d_hist_bklm_title = "Plane Efficiency Ratios in BKLM";
162  eff2d_hist_eklm_title = "Plane Efficiency Ratios in EKLM";
163  } else {
164  eff2d_hist_bklm_title = "Plane Efficiency Diffs in BKLM";
165  eff2d_hist_eklm_title = "Plane Efficiency Diffs in EKLM";
166  }
167 
168  // BKLM
169  m_eff2d_bklm = new TH2F((m_histogramDirectoryName + "/eff2d_bklm_sector").data(), eff2d_hist_bklm_title,
172  m_eff2d_bklm->GetXaxis()->SetTitle("Sector");
173  m_eff2d_bklm->GetYaxis()->SetTitle("Layer");
174  m_eff2d_bklm->SetStats(false);
175 
176  m_err_bklm = new TH2F((m_histogramDirectoryName + "/err_bklm_sector").data(), eff2d_hist_bklm_title,
179  m_err_bklm->GetXaxis()->SetTitle("Sector");
180  m_err_bklm->GetYaxis()->SetTitle("Layer");
181  m_err_bklm->SetStats(false);
182 
183  for (int sec_id = 0; sec_id < BKLMElementNumbers::getMaximalSectorNumber(); sec_id++) {
184  std::string BB_sec = "BB" + std::to_string(sec_id);
185  m_eff2d_bklm->GetXaxis()->SetBinLabel(sec_id + 1, BB_sec.c_str());
186 
187  std::string BF_sec = "BF" + std::to_string(sec_id);
188  m_eff2d_bklm->GetXaxis()->SetBinLabel(BKLMElementNumbers::getMaximalSectorNumber() + sec_id + 1, BF_sec.c_str());
189  }
190 
191  for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerNumber(); lay_id++) {
192  std::string B_lay = std::to_string(lay_id);
193  m_eff2d_bklm->GetYaxis()->SetBinLabel(lay_id + 1, B_lay.c_str());
194  }
195 
196  // EKLM
198 
199  m_eff2d_eklm = new TH2F((m_histogramDirectoryName + "/eff2d_eklm_sector").data(), eff2d_hist_eklm_title,
200  n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
202  m_eff2d_eklm->GetXaxis()->SetTitle("Sector");
203  m_eff2d_eklm->GetYaxis()->SetTitle("Layer");
204  m_eff2d_eklm->SetStats(false);
205 
206  m_err_eklm = new TH2F((m_histogramDirectoryName + "/err_eklm_sector").data(), eff2d_hist_eklm_title,
207  n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
209  m_err_eklm->GetXaxis()->SetTitle("Sector");
210  m_err_eklm->GetYaxis()->SetTitle("Layer");
211  m_err_eklm->SetStats(false);
212 
213  std::string name;
214  const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
215  for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
216  int section, layer;
218  layerGlobal, &section, &layer);
220  name = "B";
221  else
222  name = "F";
223  name += std::to_string(layer);
224  m_eff2d_eklm->GetXaxis()->SetBinLabel(layerGlobal, name.c_str());
225  }
226 
227  for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalSectorGlobalNumberKLMOrder(); lay_id++) {
228  std::string E_lay = std::to_string(lay_id);
229  m_eff2d_eklm->GetYaxis()->SetBinLabel(lay_id + 1, E_lay.c_str());
230  }
231 
232 
233 }
234 
235 
237 {
238  m_RunType = findHist("DQMInfo/rtype");
239  m_RunTypeString = m_RunType ? m_RunType->GetTitle() : "";
240  m_IsPhysicsRun = (m_RunTypeString == "physics");
241 
242  double unused = NAN;
243  //ratio/diff mode should only be possible if references exist
244  if (m_refFile && m_refFile->IsOpen()) {
245  // values for LOLO and LOW error are used for alarmThr and warnThr settings
246  // default values should be initially defined in input parameters?
247  double tempAlarm = (double) m_alarmThr;
248  double tempWarn = (double) m_warnThr;
249  requestLimitsFromEpicsPVs("2DEffSettings", tempAlarm, tempWarn, unused, unused);
250  m_alarmThr = (float) std::min(tempAlarm, tempWarn);
251  m_warnThr = (float) std::max(tempAlarm, tempWarn);
252  // EPICS should catch if this happens but just in case
253  if (m_alarmThr > m_warnThr) {
254  B2WARNING("DQMHistAnalysisKLM2Module: Found that alarmThr is greater than warnThr...");
255  }
256  }
257  m_BKLMLayerWarn = 5;
258  m_EKLMLayerWarn = 5;
259  requestLimitsFromEpicsPVs("nEffBKLMLayers", unused, unused, unused, m_BKLMLayerWarn);
260  requestLimitsFromEpicsPVs("nEffEKLMLayers", unused, unused, unused, m_EKLMLayerWarn);
261 
262 }
263 
265 {
266  std::string name;
267 
268  // Looping over the sectors
269  for (int bin = 0; bin < m_eff_bklm_sector->GetXaxis()->GetNbins(); bin++) {
270  name = "eff_B";
271  if (bin < 8)
272  name += "B";
273  else
274  name += "F";
275  name += std::to_string(bin % 8);
276  m_monObj->setVariable(name, m_eff_bklm_sector->GetBinContent(bin + 1));
277  }
278 
279  for (int bin = 0; bin < m_eff_eklm_sector->GetXaxis()->GetNbins(); bin++) {
280  name = "eff_E";
281  if (bin < 4)
282  name += "B";
283  else
284  name += "F";
285  name += std::to_string(bin % 4);
286  m_monObj->setVariable(name, m_eff_eklm_sector->GetBinContent(bin + 1));
287  }
288 
289  // Looping over the planes
290  for (int layer = 0; layer < m_eff_bklm->GetXaxis()->GetNbins(); layer++) {
291  name = "eff_B";
292  if (layer / 15 < 8) {
293  name += "B";
294  } else {
295  name += "F";
296  }
297  name += std::to_string(int(layer / 15) % 8) + "_layer" + std::to_string(1 + (layer % 15));
298  m_monObj->setVariable(name, m_eff_bklm->GetBinContent(layer + 1));
299  }
300  for (int layer = 0; layer < m_eff_eklm->GetXaxis()->GetNbins(); layer++) {
301  name = "eff_E";
302  if (layer / 8 < 12)
303  name += "B" + std::to_string(layer / 8 + 1);
304  else
305  name += "F" + std::to_string(layer / 8 - 11);
306  name += + "_num" + std::to_string(((layer) % 8) + 1);
307  m_monObj->setVariable(name, m_eff_eklm->GetBinContent(layer + 1));
308 
309  }
310 }
311 
312 void DQMHistAnalysisKLM2Module::processEfficiencyHistogram(TH1* effHist, TH1* denominator, TH1* numerator, TCanvas* canvas)
313 {
314  effHist->Reset();
315  TH1* effClone = (TH1*)effHist->Clone(); //will be useful for delta plots
316  if (denominator != nullptr && numerator != nullptr) {
317  canvas->cd();
318  effHist->Divide(numerator, denominator, 1, 1, "B");
319  effHist->Draw();
320  canvas->Modified();
321  canvas->Update();
322 
323  /* delta component */
324  auto deltaDenom = getDelta("", denominator->GetName());
325  auto deltaNumer = getDelta("", numerator->GetName());
326 
327  //both histograms should have the same update condition but checking both should be okay?
328  UpdateCanvas(canvas->GetName(), (deltaNumer != nullptr && deltaDenom != nullptr));
329  if ((deltaNumer != nullptr) && (deltaDenom != nullptr)) {
330  effClone->Divide(deltaNumer, deltaDenom, 1, 1, "B");
331  effClone->Draw("SAME");
332  canvas->Modified();
333  canvas->Update();
334  }
335  }
336 
337 }
338 
340  const std::string& histName, TH1* histogram)
341 {
342  std::string name;
343  TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
344  if (canvas == NULL) {
345  B2WARNING("KLMDQM2 histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
346  return;
347  } else {
348  canvas->Clear();
349  canvas->cd();
350  histogram->SetStats(false);
351  histogram->Draw();
352  double histMin = gPad->GetUymin();
353  double histMax = gPad->GetUymax();
354  double histRange = histMax - histMin;
355  if (histName.find("bklm") != std::string::npos) {
356  /* First draw the vertical lines and the sector names. */
357  const int maximalLayer = BKLMElementNumbers::getMaximalLayerNumber();
358  for (int sector = 0; sector < BKLMElementNumbers::getMaximalSectorGlobalNumber(); ++sector) {
359  int bin = maximalLayer * sector + 1;
360  double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
361  double xText = histogram->GetXaxis()->GetBinLowEdge(bin + maximalLayer / 2);
362  double yText = histMin + 0.98 * histRange;
363  if (sector > 0)
364  m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
365  name = "B";
366  if (sector < 8)
367  name += "B";
368  else
369  name += "F";
370  name += std::to_string(sector % 8);
371  m_PlaneText.DrawText(xText, yText, name.c_str());
372  }
373 
374  } else {
375  /* First draw the vertical lines and the sector names. */
376  const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
378  for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
379  int bin = maxPlane * layerGlobal + 1;
380  double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
381  double xText = histogram->GetXaxis()->GetBinLowEdge(bin - maxPlane / 2);
382  double yText = histMin + 0.98 * histRange;
383  if (layerGlobal < maximalLayer)
384  m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
385  int section, layer;
387  layerGlobal, &section, &layer);
389  name = "B";
390  else
391  name = "F";
392  name += std::to_string(layer);
393  m_PlaneText.DrawText(xText, yText, name.c_str());
394  }
395  }
396  canvas->Modified();
397  canvas->Update();
398  }
399 }
400 
402  TH1* mainHist, TH1* refHist, TH2* eff2dHist, TH2* errHist, int layers, int sectors, bool ratioPlot,
403  int* pvcount, double layerLimit, TCanvas* eff2dCanv)
404 {
405 
406  int i = 0;
407  float mainEff;
408  float refEff;
409  float mainErr;
410  float refErr;
411  float maxVal = m_max;
412  float minVal = m_min;
413  float eff2dVal;
414  bool setAlarm = false;
415  bool setWarn = false;
416  *pvcount = 0; //initialize to zero
417 
418  for (int binx = 0; binx < sectors; binx++) {
419 
420  for (int biny = 0; biny < layers; biny++) {
421 
422  mainEff = mainHist->GetBinContent(i + 1);
423  mainErr = mainHist->GetBinError(i + 1);
424  if (refHist) {
425  refEff = refHist->GetBinContent(i + 1);
426  refErr = refHist->GetBinError(i + 1);
427  } else {
428  refEff = 0.;
429  refErr = 0.;
430  }
431 
432  if ((mainEff == 0) and (refEff == 0)) {
433  // empty histograms, draw blank bin
434  eff2dHist->SetBinContent(binx + 1, biny + 1, 0);
435  } else if (refEff == 0) {
436  // no reference, set maximum value
437  eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
438  } else if (mainEff == 0) {
439  // no data, set zero
440  eff2dHist->SetBinContent(binx + 1, biny + 1, 0);
441  errHist->SetBinContent(binx + 1, biny + 1, 0);
442  } else {
443 
444  if (ratioPlot) {
445  eff2dVal = mainEff / refEff;
446  if (eff2dVal < m_alarmThr) {errHist->SetBinContent(binx + 1, biny + 1, eff2dVal);}
447  } else {
448  eff2dVal = (mainEff - refEff) / pow(pow(mainErr, 2) + pow(refErr, 2), 0.5);
449  }
450 
451  // main histogram
452  if ((eff2dVal > minVal) and (eff2dVal < maxVal)) {
453  // value within display window
454  eff2dHist->SetBinContent(binx + 1, biny + 1, eff2dVal);
455  } else if (eff2dVal > maxVal) {
456  // value above display window
457  eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
458  } else if (eff2dVal < minVal) {
459  // value below display window
460  eff2dHist->SetBinContent(binx + 1, biny + 1, minVal);
461  }
462 
463  // set alarm
464  if (eff2dVal < m_warnThr) {
465  *pvcount += 1;
466  }
467  if (eff2dVal < m_alarmThr) {
468  setAlarm = true;
469  }
470 
471  }
472 
473  i++;
474  }//end of layer loop
475 
476  }//end of sector loop
477 
478  if (*pvcount > (int) layerLimit) {
479  setWarn = true;
480  }
481 
482  eff2dHist->SetMinimum(m_min);
483  eff2dHist->SetMaximum(m_max);
484 
485  eff2dCanv->cd();
486  eff2dHist->Draw("COLZ");
487  errHist->Draw("TEXT SAME");
488  if (setAlarm) {
489  eff2dCanv->Pad()->SetFillColor(kRed);
490  } else if (setWarn) {
491  eff2dCanv->Pad()->SetFillColor(kYellow);
492  }
493  eff2dCanv->Modified();
494  eff2dCanv->Update();
495 }
496 
498 {
499  /* If KLM is not included, stop here and return. */
500  TH1* daqInclusion = findHist("KLM/daq_inclusion");
501  if (not(daqInclusion == NULL)) {
502  int isKlmIncluded = daqInclusion->GetBinContent(daqInclusion->GetXaxis()->FindBin("Yes"));
503  if (isKlmIncluded == 0)
504  return;
505  }
506 
507 
508  /* Obtain plots necessary for efficiency plots */
509  TH1F* all_ext_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLM");
510  TH1F* matched_hits_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLM");
511 
512  TH1F* all_ext_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLM");
513  TH1F* matched_hits_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLM");
514 
515 
516  TH1F* all_ext_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLMSector");
517  TH1F* matched_hits_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLMSector");
518 
519  TH1F* all_ext_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLMSector");
520  TH1F* matched_hits_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLMSector");
521 
522  /* Check if efficiency histograms exist*/
523  if ((all_ext_bklm == nullptr || matched_hits_bklm == nullptr) && (m_IsPhysicsRun)) {
524  B2INFO("Histograms needed for BKLM plane efficiency computation are not found");
525  }
526 
527  if ((all_ext_eklm == nullptr || matched_hits_eklm == nullptr) && (m_IsPhysicsRun)) {
528  B2INFO("Histograms needed for EKLM plane efficiency computation are not found");
529  }
530  if ((all_ext_bklm_sector == nullptr || matched_hits_bklm_sector == nullptr) && (m_IsPhysicsRun)) {
531  B2INFO("Histograms needed for BKLM sector efficiency computation are not found");
532  }
533  if ((all_ext_eklm_sector == nullptr || matched_hits_eklm_sector == nullptr) && (m_IsPhysicsRun)) {
534  B2INFO("Histograms needed for EKLM sector efficiency computation are not found");
535  }
536 
537 
538  /* Draw histograms onto canvases*/
539  processEfficiencyHistogram(m_eff_bklm, all_ext_bklm, matched_hits_bklm, m_c_eff_bklm);
540  processEfficiencyHistogram(m_eff_eklm, all_ext_eklm, matched_hits_eklm, m_c_eff_eklm);
541 
542  processEfficiencyHistogram(m_eff_bklm_sector, all_ext_bklm_sector, matched_hits_bklm_sector, m_c_eff_bklm_sector);
543  processEfficiencyHistogram(m_eff_eklm_sector, all_ext_eklm_sector, matched_hits_eklm_sector, m_c_eff_eklm_sector);
544 
545  processPlaneHistogram("eff_bklm_plane", m_eff_bklm);
546  processPlaneHistogram("eff_eklm_plane", m_eff_eklm);
547 
548  /* Make Diff 2D plots */
552 
557  /* Set EPICS PV Values*/
558  B2DEBUG(20, "Updating EPICS PVs in DQMHistAnalysisKLM2");
559  setEpicsPV("nEffBKLMLayers", m_nEffBKLMLayers);
560  setEpicsPV("nEffEKLMLayers", m_nEffEKLMLayers);
561  updateEpicsPVs(5.0);
562 }
static constexpr int getMaximalLayerGlobalNumber()
Get maximal layer global number.
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
static constexpr int getMaximalSectorNumber()
Get maximal sector number (1-based).
static constexpr int getMaximalSectorGlobalNumber()
Get maximal sector global number.
TH1 * m_ref_efficiencies_bklm
BKLM efficiencies reference histogram.
std::string m_refFileName
2D layer-sector efficiency differences
TLine m_PlaneLine
TLine for boundary in plane histograms.
int m_nEffEKLMLayers
Number of inefficient EKLM Layers.
void initialize() override final
Initializer.
TCanvas * m_c_eff2d_bklm
BKLM efficiencies ratio canvas.
double m_EKLMLayerWarn
alarm limits from inefficient EKLM layers PV
TCanvas * m_c_eff2d_eklm
EKLM efficiencies ratio canvas.
TH1 * m_eff_bklm_sector
Histogram for BKLM sector efficiency.
void processPlaneHistogram(const std::string &histName, TH1 *histogram)
Process histogram containing the number of hits in plane.
bool m_IsPhysicsRun
Run type flag for null runs.
TH1 * m_eff_eklm
Histogram for EKLM plane efficiency.
void processEfficiencyHistogram(TH1 *effHist, TH1 *denominator, TH1 *numerator, TCanvas *canvas)
Process histogram containing the efficiencies.
TCanvas * m_c_eff_eklm_sector
Histogram for EKLM sector efficiency.
const EKLMElementNumbers * m_EklmElementNumbers
EKLM element numbers.
double m_minEvents
Minimal number of entries for delta histogram update.
TH1 * m_eff_eklm_sector
Histogram for EKLM sector efficiency.
float m_max
efficiency ratio max z scale
double m_BKLMLayerWarn
alarm limits from inefficient BKLM layers PV
MonitoringObject * m_monObj
Monitoring object.
TH1 * m_ref_efficiencies_eklm
ELM efficiencies reference histogram.
float m_alarmThr
efficiency ratio alarm threshold
void process2DEffHistogram(TH1 *mainHist, TH1 *refHist, TH2 *planeHist, TH2 *errHist, int layers, int sectors, bool ratioPlot, int *pvcount, double layerLimit, TCanvas *eff2dCanv)
Process 2D efficiency histograms.
void event() override final
This method is called for each event.
TCanvas * m_c_eff_bklm_sector
Histogram for BKLM sector efficiency.
float m_warnThr
efficiency ratio warning threshold
TString m_RunTypeString
String with run type.
std::string m_histogramDirectoryName
Name of histogram directory.
TH2 * m_err_bklm
BKLM efficiencies error histogram.
float m_min
efficiency ratio min z scale
TH2 * m_eff2d_bklm
BKLM efficiencies 2dim histogram.
TH1 * m_eff_bklm
Histogram for BKLM plane efficiency.
TCanvas * m_c_eff_bklm
BKLM plane efficiency canvas.
TText m_PlaneText
TText for names in plane histograms.
TCanvas * m_c_eff_eklm
EKLM plane efficiency canvas.
void endRun() override final
This method is called if the current run ends.
void beginRun() override final
Called when entering a new run.
TH1 * m_RunType
Histogram from DQMInfo with run type.
bool m_ratio
show efficiency ratio or difference
TH2 * m_err_eklm
EKLM efficiencies error histogram.
TH2 * m_eff2d_eklm
EKLM efficiencies 2dim histogram.
int m_nEffBKLMLayers
Number of inefficient BKLM layers.
TFile * m_refFile
reference histogram file
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.
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
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)
EKLM element numbers.
static constexpr int getMaximalLayerGlobalNumber()
Get maximal detector layer global number.
void layerNumberToElementNumbers(int layerGlobal, int *section, int *layer) const
Get element numbers by detector layer global number.
static constexpr int getMaximalPlaneGlobalNumber()
Get maximal plane global number.
static constexpr int getMaximalSectorNumber()
Get maximal sector number.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
static constexpr int getMaximalSectorGlobalNumberKLMOrder()
Get maximal sector global number with KLM ordering (section, sector).
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
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.