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("KLMEfficiencyDQM"));
36 addParam("RunStopThreshold", m_stopThr, "Set stop threshold", float(0.20));
37 addParam("AlarmThreshold", m_alarmThr, "Set alarm threshold", float(0.5));
38 addParam("WarnThreshold", m_warnThr, "Set warn threshold", float(0.92));
39 addParam("Min2DEff", m_min, "2D efficiency min", float(0.5));
40 addParam("Max2DEff", m_max, "2D efficiency max", float(2));
41 addParam("RatioPlot", m_ratio, "2D efficiency ratio or difference plot ", bool(true));
42 addParam("MinEntries", m_minEntries, "Minimum entries for delta histogram update", 30000.);
43
44 m_PlaneLine.SetLineColor(kMagenta);
45 m_PlaneLine.SetLineWidth(1);
46 m_PlaneLine.SetLineStyle(2); // dashed
47 m_PlaneText.SetTextAlign(22); // centered, middle
48 m_PlaneText.SetTextColor(kMagenta);
49 m_PlaneText.SetTextFont(42); // Helvetica regular
50 m_PlaneText.SetTextSize(0.02); // 2% of TPad's full height
51}
52
53
54
56{
58
59 //register EPICS PVs
60 registerEpicsPV("KLM:Eff:nEffBKLMLayers", "nEffBKLMLayers");
61 registerEpicsPV("KLM:Eff:nEffEKLMLayers", "nEffEKLMLayers");
62 registerEpicsPV("KLM:Eff:2DEffSettings", "2DEffSettings");
63
64 gROOT->cd();
65 m_c_eff_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_plane").data());
66 m_c_eff_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_plane").data());
67 m_c_eff_bklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_bklm_sector").data());
68 m_c_eff_eklm_sector = new TCanvas((m_histogramDirectoryName + "/c_eff_eklm_sector").data());
69
70 m_c_eff2d_bklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_bklm").data());
71 m_c_eff2d_eklm = new TCanvas((m_histogramDirectoryName + "/c_eff2d_eklm").data());
72
73 m_eff_bklm = new TH1F((m_histogramDirectoryName + "/eff_bklm_plane").data(),
74 "Plane Efficiency in BKLM",
76 m_eff_bklm->GetXaxis()->SetTitle("Layer number");
77 m_eff_bklm->SetStats(false);
78 m_eff_bklm->SetOption("HIST");
79
80 m_eff_eklm = new TH1F((m_histogramDirectoryName + "/eff_eklm_plane").data(),
81 "Plane Efficiency in EKLM",
83 m_eff_eklm->GetXaxis()->SetTitle("Plane number");
84 m_eff_eklm->SetStats(false);
85 m_eff_eklm->SetOption("HIST");
86
87 m_eff_bklm_sector = new TH1F((m_histogramDirectoryName + "/eff_bklm_sector").data(),
88 "Sector Efficiency in BKLM",
90 m_eff_bklm_sector->GetXaxis()->SetTitle("Sector number");
91 m_eff_bklm_sector->SetStats(false);
92 m_eff_bklm_sector->SetOption("HIST");
93
94 m_eff_eklm_sector = new TH1F((m_histogramDirectoryName + "/eff_eklm_sector").data(),
95 "Sector Efficiency in EKLM",
97 m_eff_eklm_sector->GetXaxis()->SetTitle("Sector number");
98 m_eff_eklm_sector->SetStats(false);
99 m_eff_eklm_sector->SetOption("HIST");
100
101 /* register plots for delta histogramming */
102 // all ext hits
103 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLM", HistDelta::c_Entries, m_minEntries, 1);
104 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLM", HistDelta::c_Entries, m_minEntries, 1);
105 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsBKLMSector", HistDelta::c_Entries, m_minEntries, 1);
106 addDeltaPar(m_histogramDirectoryName, "all_ext_hitsEKLMSector", HistDelta::c_Entries, m_minEntries, 1);
107
108 // matched hits
109 addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLM", HistDelta::c_Entries, m_minEntries, 1);
110 addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLM", HistDelta::c_Entries, m_minEntries, 1);
111 addDeltaPar(m_histogramDirectoryName, "matched_hitsBKLMSector", HistDelta::c_Entries, m_minEntries, 1);
112 addDeltaPar(m_histogramDirectoryName, "matched_hitsEKLMSector", HistDelta::c_Entries, m_minEntries, 1);
113
114 // 2D Efficiency Histograms
115 TString eff2d_hist_bklm_title;
116 TString eff2d_hist_eklm_title;
117 if (m_ratio) {
118 eff2d_hist_bklm_title = "Plane Efficiency Ratios in BKLM";
119 eff2d_hist_eklm_title = "Plane Efficiency Ratios in EKLM";
120 } else {
121 eff2d_hist_bklm_title = "Plane Efficiency Diffs in BKLM";
122 eff2d_hist_eklm_title = "Plane Efficiency Diffs in EKLM";
123 }
124
125 // BKLM
126 m_eff2d_bklm = new TH2F((m_histogramDirectoryName + "/eff2d_bklm_sector").data(), eff2d_hist_bklm_title,
129 m_eff2d_bklm->GetXaxis()->SetTitle("Sector");
130 m_eff2d_bklm->GetYaxis()->SetTitle("Layer");
131 m_eff2d_bklm->SetStats(false);
132
133 m_err_bklm = new TH2F((m_histogramDirectoryName + "/err_bklm_sector").data(), eff2d_hist_bklm_title,
136 m_err_bklm->GetXaxis()->SetTitle("Sector");
137 m_err_bklm->GetYaxis()->SetTitle("Layer");
138 m_err_bklm->SetStats(false);
139
140 for (int sec_id = 0; sec_id < BKLMElementNumbers::getMaximalSectorNumber(); sec_id++) {
141 std::string BB_sec = "BB" + std::to_string(sec_id);
142 m_eff2d_bklm->GetXaxis()->SetBinLabel(sec_id + 1, BB_sec.c_str());
143
144 std::string BF_sec = "BF" + std::to_string(sec_id);
145 m_eff2d_bklm->GetXaxis()->SetBinLabel(BKLMElementNumbers::getMaximalSectorNumber() + sec_id + 1, BF_sec.c_str());
146 }
147
148 for (int lay_id = 0; lay_id < BKLMElementNumbers::getMaximalLayerNumber(); lay_id++) {
149 std::string B_lay = std::to_string(lay_id);
150 m_eff2d_bklm->GetYaxis()->SetBinLabel(lay_id + 1, B_lay.c_str());
151 }
152
153 // EKLM
155
156 m_eff2d_eklm = new TH2F((m_histogramDirectoryName + "/eff2d_eklm_sector").data(), eff2d_hist_eklm_title,
157 n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
159 m_eff2d_eklm->GetXaxis()->SetTitle("Layer");
160 m_eff2d_eklm->GetYaxis()->SetTitle("Sector");
161 m_eff2d_eklm->SetStats(false);
162
163 m_err_eklm = new TH2F((m_histogramDirectoryName + "/err_eklm_sector").data(), eff2d_hist_eklm_title,
164 n_sectors_eklm, 0.5, n_sectors_eklm + 0.5,
166 m_err_eklm->GetXaxis()->SetTitle("Layer");
167 m_err_eklm->GetYaxis()->SetTitle("Sector");
168 m_err_eklm->SetStats(false);
169
170 std::string name;
171 const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
172 for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
173 int section, layer;
175 layerGlobal, &section, &layer);
177 name = "B";
178 else
179 name = "F";
180 name += std::to_string(layer);
181 m_eff2d_eklm->GetXaxis()->SetBinLabel(layerGlobal, name.c_str());
182 }
183
184 for (int lay_id = 0; lay_id < EKLMElementNumbers::getMaximalSectorGlobalNumberKLMOrder(); lay_id++) {
185 std::string E_lay = std::to_string(lay_id);
186 m_eff2d_eklm->GetYaxis()->SetBinLabel(lay_id + 1, E_lay.c_str());
187 }
188
189
190}
191
192void DQMHistAnalysisKLM2Module::initialize2DRefHistogram(TH1*& hist, const std::string& histName, const std::string& title,
193 int maxLayer)
194{
195 // Get reference histogram using the base class method
197
198 if (hist) {
199 B2INFO("DQMHistAnalysisKLM2: " + histName + " histogram was found in reference");
200 hist->SetLineColor(2);
201 hist->SetOption("HIST");
202 hist->SetStats(false);
203 } else {
204 B2WARNING("DQMHistAnalysisKLM2: " + histName + " histogram not found in reference");
205 hist = new TH1F(histName.c_str(), title.c_str(), maxLayer, 0.5, 0.5 + maxLayer);
206 for (int lay_id = 0; lay_id < maxLayer; lay_id++) {
207 hist->SetBinContent(lay_id + 1, m_ratio ? 1 : 0);
208 }
209 }
210}
211
213{
214 m_IsPhysicsRun = (getRunType() == "physics");
215 m_IsNullRun = (getRunType() == "null");
216
217 // Initialize reference histograms for 2D efficiency plots
218 initialize2DRefHistogram(m_ref_efficiencies_bklm, "eff_bklm_plane", "Plane Efficiency in BKLM",
220
221 initialize2DRefHistogram(m_ref_efficiencies_eklm, "eff_eklm_plane", "Plane Efficiency in EKLM",
223
224 double unused = NAN;
225 //ratio/diff mode should only be possible if references exist
226 // values for LOLO, LOW & High error are used for (run-)stopThr, alarmThr and warnThr settings
227 // default values should be initially defined in input parameters?
228 double tempStop = (double) m_stopThr;
229 double tempAlarm = (double) m_alarmThr;
230 double tempWarn = (double) m_warnThr;
231 requestLimitsFromEpicsPVs("2DEffSettings", tempStop, tempAlarm, tempWarn, unused);
232
233 // Create an array of the Thresholds
234 double valuesThr[] = { tempStop, tempAlarm, tempWarn };
235
236 // Sort the array from lowest to highest
237 std::sort(std::begin(valuesThr), std::end(valuesThr));
238
239 // Assign the sorted threshold values
240 m_stopThr = (float)(valuesThr[0]); // lowest value i.e, //lolo
241 m_alarmThr = (float)(valuesThr[1]); // middle value i.e, //low
242 m_warnThr = (float)(valuesThr[2]); // highest value i.e, //high
243
244 // EPICS should catch if this happens but just in case
246 B2WARNING("DQMHistAnalysisKLM2Module: Found that alarmThr or alarmStop is greater than warnThr...");
247 }
248 m_BKLMLayerWarn = 5;
249 m_EKLMLayerWarn = 5;
250 requestLimitsFromEpicsPVs("nEffBKLMLayers", unused, unused, unused, m_BKLMLayerWarn);
251 requestLimitsFromEpicsPVs("nEffEKLMLayers", unused, unused, unused, m_EKLMLayerWarn);
252
253}
254
256{
257 std::string name;
258
259 int bklmMaxLayer = BKLMElementNumbers::getMaximalLayerNumber();//15
260 int bklmMaxSector = BKLMElementNumbers::getMaximalSectorNumber();//8
261
263 int eklmLocalMaxSector = EKLMElementNumbers::getMaximalSectorNumber();//4
265
266 // Looping over the sectors
267 for (int bin = 0; bin < m_eff_bklm_sector->GetXaxis()->GetNbins(); bin++) {
268 name = "eff_B";
269 if (bin < bklmMaxSector)
270 name += "B";
271 else
272 name += "F";
273 name += std::to_string(bin % bklmMaxSector);
274 m_monObj->setVariable(name, m_eff_bklm_sector->GetBinContent(bin + 1));
275 }
276
277 for (int bin = 0; bin < m_eff_eklm_sector->GetXaxis()->GetNbins(); bin++) {
278 name = "eff_E";
279 if (bin < eklmLocalMaxSector) //(bin < 4)
280 name += "B";
281 else
282 name += "F";
283 name += std::to_string(bin % eklmLocalMaxSector);
284 m_monObj->setVariable(name, m_eff_eklm_sector->GetBinContent(bin + 1));
285 }
286
287 // Looping over the planes
288 for (int layer = 0; layer < m_eff_bklm->GetXaxis()->GetNbins(); layer++) {
289 name = "eff_B";
290 //layer/15 < 8
291 if (layer / bklmMaxLayer < bklmMaxSector) {
292 name += "B";
293 } else {
294 name += "F";
295 }
296 name += std::to_string(int(layer / bklmMaxLayer) % bklmMaxSector) + "_layer" + std::to_string(1 + (layer % bklmMaxLayer));
297 m_monObj->setVariable(name, m_eff_bklm->GetBinContent(layer + 1));
298 }
299 for (int layer = 0; layer < m_eff_eklm->GetXaxis()->GetNbins(); layer++) {
300 name = "eff_E";
301 if (layer / eklmGlobalMaxSector < eklmBLayerCount) //(layer/8 < 12)
302 name += "B" + std::to_string(layer / eklmGlobalMaxSector + 1);
303 else
304 name += "F" + std::to_string(layer / eklmGlobalMaxSector - eklmBLayerCount + 1);
305 name += + "_num" + std::to_string(((layer) % eklmGlobalMaxSector) + 1);
306 m_monObj->setVariable(name, m_eff_eklm->GetBinContent(layer + 1));
307
308 }
309}
310
311void DQMHistAnalysisKLM2Module::processEfficiencyHistogram(TH1* effHist, TH1* denominator, TH1* numerator, TCanvas* canvas)
312{
313 effHist->Reset();
314 std::unique_ptr<TH1> effClone(static_cast<TH1*>
315 (effHist->Clone())); // Clone effHist, will be useful for delta plots & Smart pointer will manage memory leak
316 canvas->cd();
317 if (denominator != nullptr && numerator != nullptr) {
318 effHist->Divide(numerator, denominator, 1, 1, "B");
319 effHist->Draw();
320
321 //reference check
322 TH1* ref = findRefHist(effHist->GetName(), ERefScaling::c_RefScaleNone);
323 if (ref) {ref->Draw("hist,same");}
324
325 canvas->Modified();
326 canvas->Update();
327
328 /* delta component */
329 auto deltaDenom = getDelta("", denominator->GetName());
330 auto deltaNumer = getDelta("", numerator->GetName());
331
332 // both histograms should have the same update condition but checking both should be okay?
333 // if this condition is not satisfied, does it cause the above to not ever update?
334 // after test campaign, switch condition back to (deltaNumer != nullptr && deltaDenom != nullptr)
335 UpdateCanvas(canvas->GetName(), (effHist != nullptr));
336 if ((deltaNumer != nullptr) && (deltaDenom != nullptr)) {
337 B2INFO("DQMHistAnalysisKLM2: Eff Delta Num/Denom Entries is " << deltaNumer->GetEntries() << "/" << deltaDenom->GetEntries());
338 effClone->Divide(deltaNumer, deltaDenom, 1, 1, "B");
339 effClone->SetLineColor(kOrange);
340 effClone->DrawCopy("SAME"); // managed by ROOT, so it helps in plotting even if obj deleted by smart pointer
341 canvas->Modified();
342 canvas->Update();
343 }
344 }
345}
346
348 const std::string& histName, TH1* histogram)
349{
350 std::string name;
351 TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
352 if (canvas == NULL) {
353 B2WARNING("KLMDQM2 histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
354 return;
355 } else {
356 canvas->cd();
357 histogram->SetStats(false);
358 double histMin = gPad->GetUymin();
359 double histMax = gPad->GetUymax();
360 double histRange = histMax - histMin;
361 if (histName.find("bklm") != std::string::npos) {
362 /* First draw the vertical lines and the sector names. */
363 const int maximalLayer = BKLMElementNumbers::getMaximalLayerNumber();
364 for (int sector = 0; sector < BKLMElementNumbers::getMaximalSectorGlobalNumber(); ++sector) {
365 int bin = maximalLayer * sector + 1;
366 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
367 double xText = histogram->GetXaxis()->GetBinLowEdge(bin + maximalLayer / 2);
368 double yText = histMin + 0.98 * histRange;
369 if (sector > 0)
370 m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
371 name = "B";
373 name += "B";
374 else
375 name += "F";
376 name += std::to_string(sector % BKLMElementNumbers::getMaximalSectorNumber());
377 m_PlaneText.DrawText(xText, yText, name.c_str());
378 }
379
380 } else {
381 /* First draw the vertical lines and the sector names. */
382 const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
384 for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
385 int bin = maxPlane * layerGlobal + 1;
386 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
387 double xText = histogram->GetXaxis()->GetBinLowEdge(bin - maxPlane / 2);
388 double yText = histMin + 0.98 * histRange;
389 if (layerGlobal < maximalLayer)
390 m_PlaneLine.DrawLine(xLine, histMin, xLine, histMin + histRange);
391 int section, layer;
393 layerGlobal, &section, &layer);
395 name = "B";
396 else
397 name = "F";
398 name += std::to_string(layer);
399 m_PlaneText.DrawText(xText, yText, name.c_str());
400 }
401 }
402 canvas->Modified();
403 canvas->Update();
404 }
405}
406
408 TH1* mainHist, TH1* refHist, TH2* eff2dHist, TH2* errHist, int layers, int sectors, bool ratioPlot,
409 int& pvcount, double layerLimit, TCanvas* eff2dCanv)
410{
411
412 int i = 0;
413 float mainEff;
414 float refEff;
415 float mainErr;
416 float refErr;
417 float maxVal = m_max;
418 float minVal = m_min;
419 float alarmThr = m_alarmThr;
420 float warnThr = m_warnThr;
421 float stopThr = m_stopThr;
422 float eff2dVal;
423 bool setAlarm = false;
424 bool setWarn = false;
425 bool setFew = false;
426 int mainEntries;
427 int layercount = 0;
428
429 errHist->Reset(); // Reset histogram
430
431 pvcount = 0; // Initialize to zero
432 mainEntries = mainHist->GetEntries();
433
434 for (int binx = 0; binx < sectors; binx++) {
435
436 for (int biny = 0; biny < layers; biny++) {
437
438 mainEff = mainHist->GetBinContent(i + 1);
439 mainErr = mainHist->GetBinError(i + 1);
440 if (refHist) {
441 refEff = refHist->GetBinContent(i + 1);
442 refErr = refHist->GetBinError(i + 1);
443 } else {
444 refEff = 0.;
445 refErr = 0.;
446 }
447
448 if ((mainEff == 0.) and (refEff == 0.)) {
449 // empty histograms, draw blank bin
450 eff2dHist->SetBinContent(binx + 1, biny + 1, 0.);
451 } else if ((refEff == 0.) and (ratioPlot)) {
452 // no reference, set maximum value
453 eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
454 } else if (mainEff == 0.) {
455 // no data, set zero
456 eff2dHist->SetBinContent(binx + 1, biny + 1, 0.);
457 errHist->SetBinContent(binx + 1, biny + 1, 0.);
458 } else {
459
460 if (ratioPlot) {
461 eff2dVal = mainEff / refEff;
462 // Set threshold to warnThr("high") to capture efficiency drops across layers earlier
463 if (eff2dVal < warnThr) {errHist->SetBinContent(binx + 1, biny + 1, eff2dVal);}
464 } else {
465 eff2dVal = (mainEff - refEff) / pow(pow(mainErr, 2) + pow(refErr, 2), 0.5);
466 }
467
468 // main histogram
469 if ((eff2dVal > minVal) and (eff2dVal < maxVal)) {
470 // value within display window
471 eff2dHist->SetBinContent(binx + 1, biny + 1, eff2dVal);
472 } else if (eff2dVal > maxVal) {
473 // value above display window
474 eff2dHist->SetBinContent(binx + 1, biny + 1, maxVal);
475 } else if (eff2dVal < minVal) {
476 // value below display window
477 eff2dHist->SetBinContent(binx + 1, biny + 1, minVal);
478 }
479
480 // set alarm
481 if (mainEntries < (int)m_minEntries) {
482 setFew = true;
483 } else {
484 if (eff2dVal < warnThr) {
485 pvcount += 1;
486 if ((eff2dVal <= alarmThr) && (eff2dVal >= stopThr)) {
487 setWarn = true;
488 layercount++; // Increment layercount when alarm condition is met
489 } else if (eff2dVal < stopThr) {
490 setAlarm = true;
491 layercount++; // Increment layercount when a stop condition is met
492 }
493 }
494 }
495
496 }
497 i++;
498 }//end of bin y
499
500 }//end of bin x
501
502 if ((pvcount - layercount) > (int) layerLimit) {
503 setWarn = true;
504 }
505 if (layercount > (int) layerLimit) {
506 setAlarm = true;
507 }
508
509 eff2dHist->SetMinimum(minVal);
510 eff2dHist->SetMaximum(maxVal);
511
512 eff2dCanv->cd();
513 eff2dHist->Draw("COLZ");
514 errHist->Draw("TEXT SAME");
515 if (setFew) {
516 colorizeCanvas(eff2dCanv, c_StatusTooFew);
517 } else if (setAlarm) {
518 colorizeCanvas(eff2dCanv, c_StatusError);
519 } else if (setWarn) {
521 } else {
522 colorizeCanvas(eff2dCanv, c_StatusGood);
523 }
524 eff2dCanv->Modified();
525 eff2dCanv->Update();
526}
527
529{
530 /* If KLM is not included, stop here and return. */
531 TH1* daqInclusion = findHist("KLM/daq_inclusion");
532 if (not(daqInclusion == NULL)) {
533 int isKlmIncluded = daqInclusion->GetBinContent(daqInclusion->GetXaxis()->FindBin("Yes"));
534 if (isKlmIncluded == 0)
535 return;
536 }
537
538
539 /* Obtain plots necessary for efficiency plots */
540 TH1F* all_ext_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLM");
541 TH1F* matched_hits_bklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLM");
542
543 TH1F* all_ext_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLM");
544 TH1F* matched_hits_eklm = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLM");
545
546
547 TH1F* all_ext_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsBKLMSector");
548 TH1F* matched_hits_bklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsBKLMSector");
549
550 TH1F* all_ext_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/all_ext_hitsEKLMSector");
551 TH1F* matched_hits_eklm_sector = (TH1F*)findHist(m_histogramDirectoryName + "/matched_hitsEKLMSector");
552
553 /* Check if efficiency histograms exist*/
554 if ((all_ext_bklm == nullptr || matched_hits_bklm == nullptr) && (m_IsPhysicsRun)) {
555 B2INFO("Histograms needed for BKLM plane efficiency computation are not found");
556 }
557
558 if ((all_ext_eklm == nullptr || matched_hits_eklm == nullptr) && (m_IsPhysicsRun)) {
559 B2INFO("Histograms needed for EKLM plane efficiency computation are not found");
560 }
561 if ((all_ext_bklm_sector == nullptr || matched_hits_bklm_sector == nullptr) && (m_IsPhysicsRun)) {
562 B2INFO("Histograms needed for BKLM sector efficiency computation are not found");
563 }
564 if ((all_ext_eklm_sector == nullptr || matched_hits_eklm_sector == nullptr) && (m_IsPhysicsRun)) {
565 B2INFO("Histograms needed for EKLM sector efficiency computation are not found");
566 }
567
568
569 /* Draw histograms onto canvases*/
570 processEfficiencyHistogram(m_eff_bklm, all_ext_bklm, matched_hits_bklm, m_c_eff_bklm);
571 processEfficiencyHistogram(m_eff_eklm, all_ext_eklm, matched_hits_eklm, m_c_eff_eklm);
572
573 processEfficiencyHistogram(m_eff_bklm_sector, all_ext_bklm_sector, matched_hits_bklm_sector, m_c_eff_bklm_sector);
574 processEfficiencyHistogram(m_eff_eklm_sector, all_ext_eklm_sector, matched_hits_eklm_sector, m_c_eff_eklm_sector);
575
576 processPlaneHistogram("eff_bklm_plane", m_eff_bklm);
577 processPlaneHistogram("eff_eklm_plane", m_eff_eklm);
578
579 /* Make Diff 2D plots */
583
588 /* Set EPICS PV Values*/
589 B2DEBUG(20, "DQMHistAnalysisKLM2: Updating EPICS PVs");
590 // only update PVs if there's enough statistics and datasize != 0
591 // Check if it's a null run, if so, don't update EPICS PVs
592 if (m_IsNullRun) {
593 B2INFO("DQMHistAnalysisKLM2: Null run detected. No PV Update.");
594 return;
595 }
596 auto* daqDataSize = findHist("DAQ/KLMDataSize");
597 double meanDAQDataSize = 0.;
598 if (daqDataSize != nullptr) {
599 meanDAQDataSize = daqDataSize->GetMean();
600 B2INFO("DAQ/KLMDataSize's mean is " << meanDAQDataSize);
601 } else
602 B2WARNING("DQMHistAnalysisKLM2: Cannot find KLMDataSize");
603 if ((daqDataSize != nullptr) and (meanDAQDataSize != 0.)) {
604 int procesedEvents = DQMHistAnalysisModule::getEventProcessed();
605 if (procesedEvents > (int)m_minEntries) {
606 if (static_cast<int>(m_eff_bklm->GetEntries()) > (int)m_minEntries) {
607 setEpicsPV("nEffBKLMLayers", m_nEffBKLMLayers);
608 }
609 if (static_cast<int>(m_eff_eklm->GetEntries()) > (int)m_minEntries) {
610 setEpicsPV("nEffEKLMLayers", m_nEffEKLMLayers);
611 }
612 }
613 } else
614 B2INFO("DQMHistAnalysisKLM2: KLM Not included. No PV Update. ");
615
616}
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.
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
warn 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
warn 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 event() override final
This method is called for each event.
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.
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.
void initialize2DRefHistogram(TH1 *&hist, const std::string &histName, const std::string &title, int maxLayer)
Initialize a histogram by either finding a reference histogram or creating a new one.
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.
The base class for the histogram analysis module.
TCanvas * findCanvas(TString cname)
Find canvas by name.
static TH1 * findRefHist(const std::string &histname, ERefScaling scaling=ERefScaling::c_RefScaleNone, const TH1 *hist=nullptr)
Get referencehistogram from list (no other search).
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.
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.