Belle II Software development
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
23using namespace Belle2;
24
25REG_MODULE(DQMHistAnalysisKLM2);
26
29 m_IsNullRun{false},
30 m_EklmElementNumbers{&(EKLMElementNumbers::Instance())}
31{
32 setDescription("Module used to analyze KLM Efficiency DQM histograms (depends on tracking variables).");
33 addParam("HistogramDirectoryName", m_histogramDirectoryName, "Name of histogram directory", std::string("KLMEfficiencyDQM"));
34 addParam("RefHistogramDirectoryName", m_refHistogramDirectoryName, "Name of ref histogram directory",
35 std::string("ref/KLMEfficiencyDQM"));
36 addParam("RefHistoFile", m_refFileName, "Reference histogram file name", std::string("KLM_DQM_REF_BEAM.root"));
37 addParam("RunStopThreshold", m_stopThr, "Set stop threshold", float(0.20));
38 addParam("AlarmThreshold", m_alarmThr, "Set alarm threshold", float(0.5));
39 addParam("WarnThreshold", m_warnThr, "Set warn threshold", float(0.92));
40 addParam("Min2DEff", m_min, "2D efficiency min", float(0.5));
41 addParam("Max2DEff", m_max, "2D efficiency max", float(2));
42 addParam("RatioPlot", m_ratio, "2D efficiency ratio or difference plot ", bool(true));
43 addParam("MinEntries", m_minEntries, "Minimum entries for delta histogram update", 30000.);
44
45 m_PlaneLine.SetLineColor(kMagenta);
46 m_PlaneLine.SetLineWidth(1);
47 m_PlaneLine.SetLineStyle(2); // dashed
48 m_PlaneText.SetTextAlign(22); // centered, middle
49 m_PlaneText.SetTextColor(kMagenta);
50 m_PlaneText.SetTextFont(42); // Helvetica regular
51 m_PlaneText.SetTextSize(0.02); // 2% of TPad's full height
52}
53
54
55
57{
59
60 //register EPICS PVs
61 registerEpicsPV("KLM:Eff:nEffBKLMLayers", "nEffBKLMLayers");
62 registerEpicsPV("KLM:Eff:nEffEKLMLayers", "nEffEKLMLayers");
63 registerEpicsPV("KLM:Eff:2DEffSettings", "2DEffSettings");
64
65 if (m_refFileName != "") {
66 m_refFile = TFile::Open(m_refFileName.data(), "READ");
67 }
68
69 // Get Reference Histograms
70 if (m_refFile && m_refFile->IsOpen()) {
71 B2INFO("DQMHistAnalysisKLM2: reference root file (" << m_refFileName << ") FOUND, able to read ref histograms");
72
73 m_ref_efficiencies_bklm = (TH1F*)m_refFile->Get((m_refHistogramDirectoryName + "/eff_bklm_plane").data());
74 if (m_ref_efficiencies_bklm != nullptr) {
75 B2INFO("DQMHistAnalysisKLM2: eff_bklm_plane histogram was found in reference");
76 m_ref_efficiencies_bklm->SetLineColor(2);
77 m_ref_efficiencies_bklm->SetOption("HIST");
78 m_ref_efficiencies_bklm->SetStats(false);
79 } else {
80 B2WARNING("DQMHistAnalysisKLM2: eff_bklm_plane histogram not found in reference");
81 m_ref_efficiencies_bklm = new TH1F("eff_bklm_plane", "Plane Efficiency in BKLM", BKLMElementNumbers::getMaximalLayerGlobalNumber(),
83 for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerGlobalNumber(); lay_id++) {
84 if (m_ratio) {
85 m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 1);
86 } else {
87 m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 0);
88 }
89 }
90 }
91
92
93 m_ref_efficiencies_eklm = (TH1F*)m_refFile->Get((m_refHistogramDirectoryName + "/eff_eklm_plane").data());
94 if (m_ref_efficiencies_eklm != nullptr) {
95 B2INFO("DQMHistAnalysisKLM2: eff_eklm_plane histogram was found in reference");
96 m_ref_efficiencies_eklm->SetLineColor(2);
97 m_ref_efficiencies_eklm->SetOption("HIST");
98 m_ref_efficiencies_eklm->SetStats(false);
99 } else {
100 B2WARNING("DQMHistAnalysisKLM2: eff_eklm_plane histogram not found in reference");
101 m_ref_efficiencies_eklm = new TH1F("eff_eklm_plane", "Plane Efficiency in EKLM", EKLMElementNumbers::getMaximalPlaneGlobalNumber(),
103 for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalPlaneGlobalNumber(); lay_id++) {
104 if (m_ratio) {
105 m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 1);
106 } else {
107 m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 0);
108 }
109 }
110 }
111 } else {
112 B2WARNING("DQMHistAnalysisKLM2: reference root file (" << m_refFileName << ") not found, or closed");
113
114 // Switch to absolute 2D efficiencies if reference histogram is not found
115 m_stopThr = 0.0;
116 m_alarmThr = 0.35;
117 m_warnThr = 0.5; //contigency value to still spot some problems
118 m_ref_efficiencies_bklm = new TH1F("eff_bklm_plane", "Plane Efficiency in BKLM", BKLMElementNumbers::getMaximalLayerGlobalNumber(),
120 for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerGlobalNumber(); lay_id++) {
121 if (m_ratio) {
122 m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 1);
123 } else {
124 m_ref_efficiencies_bklm->SetBinContent(lay_id + 1, 0);
125 }
126 }
127
128 m_ref_efficiencies_eklm = new TH1F("eff_eklm_plane", "Plane Efficiency in EKLM", EKLMElementNumbers::getMaximalPlaneGlobalNumber(),
130 for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalPlaneGlobalNumber(); lay_id++) {
131 if (m_ratio) {
132 m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 1);
133 } else {
134 m_ref_efficiencies_eklm->SetBinContent(lay_id + 1, 0);
135 }
136 }
137
138 }
139 gROOT->cd();
140 m_c_eff_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_plane").data());
141 m_c_eff_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_plane").data());
142 m_c_eff_bklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_sector").data());
143 m_c_eff_eklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_sector").data());
144
145 m_c_eff2d_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_bklm").data());
146 m_c_eff2d_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_eklm").data());
147
148 m_eff_bklm = new TH1F((m_histogramDirectoryName + "/eff_bklm_plane").data(),
149 "Plane Efficiency in BKLM",
151 m_eff_bklm->GetXaxis()->SetTitle("Layer number");
152 m_eff_bklm->SetStats(false);
153 m_eff_bklm->SetOption("HIST");
154
155 m_eff_eklm = new TH1F((m_histogramDirectoryName + "/eff_eklm_plane").data(),
156 "Plane Efficiency in EKLM",
158 m_eff_eklm->GetXaxis()->SetTitle("Plane number");
159 m_eff_eklm->SetStats(false);
160 m_eff_eklm->SetOption("HIST");
161
162 m_eff_bklm_sector = new TH1F((m_histogramDirectoryName + "/eff_bklm_sector").data(),
163 "Sector Efficiency in BKLM",
165 m_eff_bklm_sector->GetXaxis()->SetTitle("Sector number");
166 m_eff_bklm_sector->SetStats(false);
167 m_eff_bklm_sector->SetOption("HIST");
168
169 m_eff_eklm_sector = new TH1F((m_histogramDirectoryName + "/eff_eklm_sector").data(),
170 "Sector Efficiency in EKLM",
172 m_eff_eklm_sector->GetXaxis()->SetTitle("Sector number");
173 m_eff_eklm_sector->SetStats(false);
174 m_eff_eklm_sector->SetOption("HIST");
175
176 /* register plots for delta histogramming */
177 // all ext hits
178 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLM", HistDelta::c_Entries, m_minEntries, 1);
179 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLM", HistDelta::c_Entries, m_minEntries, 1);
180 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLMSector", HistDelta::c_Entries, m_minEntries, 1);
181 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLMSector", HistDelta::c_Entries, m_minEntries, 1);
182
183 // matched hits
184 addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLM", HistDelta::c_Entries, m_minEntries, 1);
185 addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLM", HistDelta::c_Entries, m_minEntries, 1);
186 addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLMSector", HistDelta::c_Entries, m_minEntries, 1);
187 addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLMSector", HistDelta::c_Entries, m_minEntries, 1);
188
189 // 2D Efficiency Histograms
190 TString eff2d_hist_bklm_title;
191 TString eff2d_hist_eklm_title;
192 if (m_ratio) {
193 eff2d_hist_bklm_title = "Plane Efficiency Ratios in BKLM";
194 eff2d_hist_eklm_title = "Plane Efficiency Ratios in EKLM";
195 } else {
196 eff2d_hist_bklm_title = "Plane Efficiency Diffs in BKLM";
197 eff2d_hist_eklm_title = "Plane Efficiency Diffs in EKLM";
198 }
199
200 // BKLM
201 m_eff2d_bklm = new TH2F((m_histogramDirectoryName + "/eff2d_bklm_sector").data(), eff2d_hist_bklm_title,
204 m_eff2d_bklm->GetXaxis()->SetTitle("Sector");
205 m_eff2d_bklm->GetYaxis()->SetTitle("Layer");
206 m_eff2d_bklm->SetStats(false);
207
208 m_err_bklm = new TH2F((m_histogramDirectoryName + "/err_bklm_sector").data(), eff2d_hist_bklm_title,
211 m_err_bklm->GetXaxis()->SetTitle("Sector");
212 m_err_bklm->GetYaxis()->SetTitle("Layer");
213 m_err_bklm->SetStats(false);
214
215 for (int sec_id = 0; sec_id < BKLMElementNumbers::getMaximalSectorNumber(); sec_id++) {
216 std::string BB_sec = "BB" + std::to_string(sec_id);
217 m_eff2d_bklm->GetXaxis()->SetBinLabel(sec_id + 1, BB_sec.c_str());
218
219 std::string BF_sec = "BF" + std::to_string(sec_id);
220 m_eff2d_bklm->GetXaxis()->SetBinLabel(BKLMElementNumbers::getMaximalSectorNumber() + sec_id + 1, BF_sec.c_str());
221 }
222
223 for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerNumber(); lay_id++) {
224 std::string B_lay = std::to_string(lay_id);
225 m_eff2d_bklm->GetYaxis()->SetBinLabel(lay_id + 1, B_lay.c_str());
226 }
227
228 // EKLM
230
231 m_eff2d_eklm = new TH2F((m_histogramDirectoryName + "/eff2d_eklm_sector").data(), eff2d_hist_eklm_title,
232 n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
234 m_eff2d_eklm->GetXaxis()->SetTitle("Layer");
235 m_eff2d_eklm->GetYaxis()->SetTitle("Sector");
236 m_eff2d_eklm->SetStats(false);
237
238 m_err_eklm = new TH2F((m_histogramDirectoryName + "/err_eklm_sector").data(), eff2d_hist_eklm_title,
239 n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
241 m_err_eklm->GetXaxis()->SetTitle("Layer");
242 m_err_eklm->GetYaxis()->SetTitle("Sector");
243 m_err_eklm->SetStats(false);
244
245 std::string name;
246 const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
247 for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
248 int section, layer;
250 layerGlobal, &section, &layer);
252 name = "B";
253 else
254 name = "F";
255 name += std::to_string(layer);
256 m_eff2d_eklm->GetXaxis()->SetBinLabel(layerGlobal, name.c_str());
257 }
258
259 for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalSectorGlobalNumberKLMOrder(); lay_id++) {
260 std::string E_lay = std::to_string(lay_id);
261 m_eff2d_eklm->GetYaxis()->SetBinLabel(lay_id + 1, E_lay.c_str());
262 }
263
264
265}
266
267
269{
270 m_IsPhysicsRun = (getRunType() == "physics");
271 m_IsNullRun = (getRunType() == "null");
272
273 double unused = NAN;
274 //ratio/diff mode should only be possible if references exist
275 if (m_refFile && m_refFile->IsOpen()) {
276 // values for LOLO, LOW & High error are used for (run-)stopThr, alarmThr and warnThr settings
277 // default values should be initially defined in input parameters?
278 double tempStop = (double) m_stopThr;
279 double tempAlarm = (double) m_alarmThr;
280 double tempWarn = (double) m_warnThr;
281 requestLimitsFromEpicsPVs("2DEffSettings", tempStop, tempAlarm, tempWarn, unused);
282
283 // Create an array of the Thresholds
284 double valuesThr[] = { tempStop, tempAlarm, tempWarn };
285
286 // Sort the array from lowest to highest
287 std::sort(std::begin(valuesThr), std::end(valuesThr));
288
289 // Assign the sorted threshold values
290 m_stopThr = (float)(valuesThr[0]); // lowest value i.e, //lolo
291 m_alarmThr = (float)(valuesThr[1]); // middle value i.e, //low
292 m_warnThr = (float)(valuesThr[2]); // highest value i.e, //high
293
294 // EPICS should catch if this happens but just in case
296 B2WARNING("DQMHistAnalysisKLM2Module: Found that alarmThr or alarmStop is greater than warnThr...");
297 }
298 }
299 m_BKLMLayerWarn = 5;
300 m_EKLMLayerWarn = 5;
301 requestLimitsFromEpicsPVs("nEffBKLMLayers", unused, unused, unused, m_BKLMLayerWarn);
302 requestLimitsFromEpicsPVs("nEffEKLMLayers", unused, unused, unused, m_EKLMLayerWarn);
303
304}
305
307{
308 std::string name;
309
310 int bklmMaxLayer = BKLMElementNumbers::getMaximalLayerNumber();//15
311 int bklmMaxSector = BKLMElementNumbers::getMaximalSectorNumber();//8
312
314 int eklmLocalMaxSector = EKLMElementNumbers::getMaximalSectorNumber();//4
316
317 // Looping over the sectors
318 for (int bin = 0; bin < m_eff_bklm_sector->GetXaxis()->GetNbins(); bin++) {
319 name = "eff_B";
320 if (bin < bklmMaxSector)
321 name += "B";
322 else
323 name += "F";
324 name += std::to_string(bin % bklmMaxSector);
325 m_monObj->setVariable(name, m_eff_bklm_sector->GetBinContent(bin + 1));
326 }
327
328 for (int bin = 0; bin < m_eff_eklm_sector->GetXaxis()->GetNbins(); bin++) {
329 name = "eff_E";
330 if (bin < eklmLocalMaxSector) //(bin < 4)
331 name += "B";
332 else
333 name += "F";
334 name += std::to_string(bin % eklmLocalMaxSector);
335 m_monObj->setVariable(name, m_eff_eklm_sector->GetBinContent(bin + 1));
336 }
337
338 // Looping over the planes
339 for (int layer = 0; layer < m_eff_bklm->GetXaxis()->GetNbins(); layer++) {
340 name = "eff_B";
341 //layer/15 < 8
342 if (layer / bklmMaxLayer < bklmMaxSector) {
343 name += "B";
344 } else {
345 name += "F";
346 }
347 name += std::to_string(int(layer / bklmMaxLayer) % bklmMaxSector) + "_layer" + std::to_string(1 + (layer % bklmMaxLayer));
348 m_monObj->setVariable(name, m_eff_bklm->GetBinContent(layer + 1));
349 }
350 for (int layer = 0; layer < m_eff_eklm->GetXaxis()->GetNbins(); layer++) {
351 name = "eff_E";
352 if (layer / eklmGlobalMaxSector < eklmBLayerCount) //(layer/8 < 12)
353 name += "B" + std::to_string(layer / eklmGlobalMaxSector + 1);
354 else
355 name += "F" + std::to_string(layer / eklmGlobalMaxSector - eklmBLayerCount + 1);
356 name += + "_num" + std::to_string(((layer) % eklmGlobalMaxSector) + 1);
357 m_monObj->setVariable(name, m_eff_eklm->GetBinContent(layer + 1));
358
359 }
360}
361
362void DQMHistAnalysisKLM2Module::processEfficiencyHistogram(TH1* effHist, TH1* denominator, TH1* numerator, TCanvas* canvas)
363{
364 effHist->Reset();
365 std::unique_ptr<TH1> effClone(static_cast<TH1*>
366 (effHist->Clone())); // Clone effHist, will be useful for delta plots & Smart pointer will manage memory leak
367 canvas->cd();
368 if (denominator != nullptr && numerator != nullptr) {
369 effHist->Divide(numerator, denominator, 1, 1, "B");
370 effHist->Draw();
371
372 //reference check
373 TH1* ref = findRefHist(effHist->GetName(), false);
374 if (ref) {ref->Draw("hist,same");}
375
376 canvas->Modified();
377 canvas->Update();
378
379 /* delta component */
380 auto deltaDenom = getDelta("", denominator->GetName());
381 auto deltaNumer = getDelta("", numerator->GetName());
382
383 // both histograms should have the same update condition but checking both should be okay?
384 // if this condition is not satisfied, does it cause the above to not ever update?
385 // after test campaign, switch condition back to (deltaNumer != nullptr && deltaDenom != nullptr)
386 UpdateCanvas(canvas->GetName(), (effHist != nullptr));
387 if ((deltaNumer != nullptr) && (deltaDenom != nullptr)) {
388 B2INFO("DQMHistAnalysisKLM2: Eff Delta Num/Denom Entries is " << deltaNumer->GetEntries() << "/" << deltaDenom->GetEntries());
389 effClone->Divide(deltaNumer, deltaDenom, 1, 1, "B");
390 effClone->SetLineColor(kOrange);
391 effClone->DrawCopy("SAME"); // managed by ROOT, so it helps in plotting even if obj deleted by smart pointer
392 canvas->Modified();
393 canvas->Update();
394 }
395 }
396}
397
399 const std::string& histName, TH1* histogram)
400{
401 std::string name;
402 TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
403 if (canvas == NULL) {
404 B2WARNING("KLMDQM2 histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
405 return;
406 } else {
407 canvas->cd();
408 histogram->SetStats(false);
409 double histMin = gPad->GetUymin();
410 double histMax = gPad->GetUymax();
411 double histRange = histMax - histMin;
412 if (histName.find("bklm") != std::string::npos) {
413 /* First draw the vertical lines and the sector names. */
414 const int maximalLayer = BKLMElementNumbers::getMaximalLayerNumber();
415 for (int sector = 0; sector < BKLMElementNumbers::getMaximalSectorGlobalNumber(); ++sector) {
416 int bin = maximalLayer * sector + 1;
417 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
418 double xText = histogram->GetXaxis()->GetBinLowEdge(bin + maximalLayer / 2);
419 double yText = histMin + 0.98 * histRange;
420 if (sector > 0)
421 m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
422 name = "B";
424 name += "B";
425 else
426 name += "F";
427 name += std::to_string(sector % BKLMElementNumbers::getMaximalSectorNumber());
428 m_PlaneText.DrawText(xText, yText, name.c_str());
429 }
430
431 } else {
432 /* First draw the vertical lines and the sector names. */
433 const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
435 for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
436 int bin = maxPlane * layerGlobal + 1;
437 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
438 double xText = histogram->GetXaxis()->GetBinLowEdge(bin - maxPlane / 2);
439 double yText = histMin + 0.98 * histRange;
440 if (layerGlobal < maximalLayer)
441 m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
442 int section, layer;
444 layerGlobal, &section, &layer);
446 name = "B";
447 else
448 name = "F";
449 name += std::to_string(layer);
450 m_PlaneText.DrawText(xText, yText, name.c_str());
451 }
452 }
453 canvas->Modified();
454 canvas->Update();
455 }
456}
457
459 TH1* mainHist, TH1* refHist, TH2* eff2dHist, TH2* errHist, int layers, int sectors, bool ratioPlot,
460 int* pvcount, double layerLimit, TCanvas* eff2dCanv)
461{
462
463 int i = 0;
464 float mainEff;
465 float refEff;
466 float mainErr;
467 float refErr;
468 float maxVal = m_max;
469 float minVal = m_min;
470 float alarmThr = m_alarmThr;
471 float warnThr = m_warnThr;
472 float stopThr = m_stopThr;
473 float eff2dVal;
474 bool setAlarm = false;
475 bool setWarn = false;
476 bool setFew = false;
477 int mainEntries;
478
479 errHist->Reset(); // Reset histogram
480
481 *pvcount = 0; //initialize to zero
482 mainEntries = mainHist->GetEntries();
483
484 for (int binx = 0; binx < sectors; binx++) {
485
486 for (int biny = 0; biny < layers; biny++) {
487
488 mainEff = mainHist->GetBinContent(i + 1);
489 mainErr = mainHist->GetBinError(i + 1);
490 if (refHist) {
491 refEff = refHist->GetBinContent(i + 1);
492 refErr = refHist->GetBinError(i + 1);
493 } else {
494 refEff = 0.;
495 refErr = 0.;
496 }
497
498 if ((mainEff == 0.) and (refEff == 0.)) {
499 // empty histograms, draw blank bin
500 eff2dHist->SetBinContent(binx + 1, biny + 1, 0.);
501 } else if ((refEff == 0.) and (ratioPlot)) {
502 // no reference, set maximum value
503 eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
504 } else if (mainEff == 0.) {
505 // no data, set zero
506 eff2dHist->SetBinContent(binx + 1, biny + 1, 0.);
507 errHist->SetBinContent(binx + 1, biny + 1, 0.);
508 } else {
509
510 if (ratioPlot) {
511 eff2dVal = mainEff / refEff;
512 if (eff2dVal < alarmThr) {errHist->SetBinContent(binx + 1, biny + 1, eff2dVal);}
513 } else {
514 eff2dVal = (mainEff - refEff) / pow(pow(mainErr, 2) + pow(refErr, 2), 0.5);
515 }
516
517 // main histogram
518 if ((eff2dVal > minVal) and (eff2dVal < maxVal)) {
519 // value within display window
520 eff2dHist->SetBinContent(binx + 1, biny + 1, eff2dVal);
521 } else if (eff2dVal > maxVal) {
522 // value above display window
523 eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
524 } else if (eff2dVal < minVal) {
525 // value below display window
526 eff2dHist->SetBinContent(binx + 1, biny + 1, minVal);
527 }
528
529 // set alarm
530 if (mainEntries < (int)m_minEntries) {
531 setFew = true;
532 } else {
533 if (eff2dVal < warnThr) {
534 *pvcount += 1;
535 if ((eff2dVal <= alarmThr) && (eff2dVal >= stopThr)) {
536 setWarn = true;
537 } else if (eff2dVal < stopThr) {
538 setAlarm = true;
539 }
540 }
541 }
542
543 }
544 i++;
545 }//end of bin y
546
547 }//end of bin x
548
549 if (*pvcount > (int) layerLimit) {
550 setAlarm = true;
551 }
552
553 eff2dHist->SetMinimum(minVal);
554 eff2dHist->SetMaximum(maxVal);
555
556 eff2dCanv->cd();
557 eff2dHist->Draw("COLZ");
558 errHist->Draw("TEXT SAME");
559 if (setFew) {
560 colorizeCanvas(eff2dCanv, c_StatusTooFew);
561 } else if (setAlarm) {
562 colorizeCanvas(eff2dCanv, c_StatusError);
563 } else if (setWarn) {
565 } else {
566 colorizeCanvas(eff2dCanv, c_StatusGood);
567 }
568 eff2dCanv->Modified();
569 eff2dCanv->Update();
570}
571
573{
574 /* If KLM is not included, stop here and return. */
575 TH1* daqInclusion = findHist("KLM/daq_inclusion");
576 if (not(daqInclusion == NULL)) {
577 int isKlmIncluded = daqInclusion->GetBinContent(daqInclusion->GetXaxis()->FindBin("Yes"));
578 if (isKlmIncluded == 0)
579 return;
580 }
581
582
583 /* Obtain plots necessary for efficiency plots */
584 TH1F* all_ext_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLM");
585 TH1F* matched_hits_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLM");
586
587 TH1F* all_ext_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLM");
588 TH1F* matched_hits_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLM");
589
590
591 TH1F* all_ext_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLMSector");
592 TH1F* matched_hits_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLMSector");
593
594 TH1F* all_ext_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLMSector");
595 TH1F* matched_hits_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLMSector");
596
597 /* Check if efficiency histograms exist*/
598 if ((all_ext_bklm == nullptr || matched_hits_bklm == nullptr) && (m_IsPhysicsRun)) {
599 B2INFO("Histograms needed for BKLM plane efficiency computation are not found");
600 }
601
602 if ((all_ext_eklm == nullptr || matched_hits_eklm == nullptr) && (m_IsPhysicsRun)) {
603 B2INFO("Histograms needed for EKLM plane efficiency computation are not found");
604 }
605 if ((all_ext_bklm_sector == nullptr || matched_hits_bklm_sector == nullptr) && (m_IsPhysicsRun)) {
606 B2INFO("Histograms needed for BKLM sector efficiency computation are not found");
607 }
608 if ((all_ext_eklm_sector == nullptr || matched_hits_eklm_sector == nullptr) && (m_IsPhysicsRun)) {
609 B2INFO("Histograms needed for EKLM sector efficiency computation are not found");
610 }
611
612
613 /* Draw histograms onto canvases*/
614 processEfficiencyHistogram(m_eff_bklm, all_ext_bklm, matched_hits_bklm, m_c_eff_bklm);
615 processEfficiencyHistogram(m_eff_eklm, all_ext_eklm, matched_hits_eklm, m_c_eff_eklm);
616
617 processEfficiencyHistogram(m_eff_bklm_sector, all_ext_bklm_sector, matched_hits_bklm_sector, m_c_eff_bklm_sector);
618 processEfficiencyHistogram(m_eff_eklm_sector, all_ext_eklm_sector, matched_hits_eklm_sector, m_c_eff_eklm_sector);
619
620 processPlaneHistogram("eff_bklm_plane", m_eff_bklm);
621 processPlaneHistogram("eff_eklm_plane", m_eff_eklm);
622
623 /* Make Diff 2D plots */
627
632 /* Set EPICS PV Values*/
633 B2DEBUG(20, "DQMHistAnalysisKLM2: Updating EPICS PVs");
634 // only update PVs if there's enough statistics and datasize != 0
635 // Check if it's a null run, if so, don't update EPICS PVs
636 if (m_IsNullRun) {
637 B2INFO("DQMHistAnalysisKLM2: Null run detected. No PV Update.");
638 return;
639 }
640 auto* daqDataSize = findHist("DAQ/KLMDataSize");
641 double meanDAQDataSize = 0.;
642 if (daqDataSize != nullptr) {
643 meanDAQDataSize = daqDataSize->GetMean();
644 B2INFO("DAQ/KLMDataSize's mean is " << meanDAQDataSize);
645 } else
646 B2WARNING("DQMHistAnalysisKLM2: Cannot find KLMDataSize");
647 if ((daqDataSize != nullptr) and (meanDAQDataSize != 0.)) {
648 int procesedEvents = DQMHistAnalysisModule::getEventProcessed();
649 if (procesedEvents > (int)m_minEntries) {
650 if (static_cast<int>(m_eff_bklm->GetEntries()) > (int)m_minEntries) {
651 setEpicsPV("nEffBKLMLayers", m_nEffBKLMLayers);
652 }
653 if (static_cast<int>(m_eff_eklm->GetEntries()) > (int)m_minEntries) {
654 setEpicsPV("nEffEKLMLayers", m_nEffEKLMLayers);
655 }
656 }
657 } else
658 B2INFO("DQMHistAnalysisKLM2: KLM Not included. No PV Update. ");
659
660}
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
std::string m_refHistogramDirectoryName
Name of histogram directory for reference file.
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.
float m_stopThr
efficiency ratio (run-)stop threshold
double m_minEntries
Minimal number of entries for delta histogram and PV update.
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.
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
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.
bool m_ratio
show efficiency ratio or difference
TH2 * m_err_eklm
EKLM efficiencies error histogram.
bool m_IsNullRun
Run type flag for null runs.
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.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
void colorizeCanvas(TCanvas *canvas, EStatus status)
Helper function for Canvas colorization.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static const std::string & getRunType(void)
Get the Run Type.
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
static TH1 * findRefHist(const std::string &histname, int scaling=0, const TH1 *hist=nullptr)
Get referencehistogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
@ c_StatusTooFew
Not enough entries/event to judge.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusWarning
Analysis result: Warning, there may be minor issues.
@ c_StatusGood
Analysis result: Good.
static int getEventProcessed(void)
Get the number of processed events.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
EKLM element numbers.
static constexpr int getMaximalLayerGlobalNumber()
Get maximal detector layer global number.
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer 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)
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
#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.