Belle II Software development
DQMHistAnalysisPXDEff.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// File : DQMHistAnalysisPXDEff.cc
10// Description : DQM module, which gives histograms showing the efficiency of PXD sensors
11//-
12
13#include <dqm/analysis/modules/DQMHistAnalysisPXDEff.h>
14#include <TROOT.h>
15#include <TLatex.h>
16#include <TGraphAsymmErrors.h>
17#include <vxd/geometry/GeoCache.h>
18
19using namespace std;
20using namespace Belle2;
21
22//-----------------------------------------------------------------
23// Register the Module
24//-----------------------------------------------------------------
25REG_MODULE(DQMHistAnalysisPXDEff);
26
27//-----------------------------------------------------------------
28// Implementation
29//-----------------------------------------------------------------
30
32{
33 // This module CAN NOT be run in parallel!
34 setDescription("DQM Analysis for PXD Efficiency");
35
36 // Parameter definition
37
38 // Would be much more elegant to get bin numbers from the saved histograms, but would need to retrieve at least one of them before the initialize function for this
39 // Or get one and clone it
40 addParam("binsU", m_u_bins, "histogram bins in u direction, needs to be the same as in PXDDQMEfficiency", int(16));
41 addParam("binsV", m_v_bins, "histogram bins in v direction, needs to be the same as in PXDDQMEfficiency", int(48));
42 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms were placed",
43 std::string("PXDEFF"));
44 addParam("ConfidenceLevel", m_confidence, "Confidence Level for error bars and alarms", 0.99);
45 addParam("WarnLevel", m_warnlevel, "Efficiency Warn Level for alarms", 0.92);
46 addParam("ErrorLevel", m_errorlevel, "Efficiency Level for alarms", 0.90);
47 addParam("perModuleAlarm", m_perModuleAlarm, "Alarm level per module", true);
48 addParam("alarmAdhoc", m_alarmAdhoc, "Generate Alarm from adhoc values", true);
49 addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 1000);
50 addParam("excluded", m_excluded, "the list of excluded modules, indices from 0 to 39");
51 B2DEBUG(1, "DQMHistAnalysisPXDEff: Constructor done.");
52}
53
55{
56 B2DEBUG(99, "DQMHistAnalysisPXDEffModule: initialized.");
57
60
61 // collect the list of all PXD Modules in the geometry here
62 std::vector<VxdID> sensors = geo.getListOfSensors();
63 for (VxdID& aVxdID : sensors) {
64 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
65 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
66 m_PXDModules.push_back(aVxdID); // reorder, sort would be better
67 }
68 std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
69
70 gROOT->cd(); // this seems to be important, or strange things happen
71
72 int nu = 1;//If this does not get overwritten, the histograms will anyway never contain anything useful
73 int nv = 1;
74 if (m_PXDModules.size() == 0) {
75 // This could as well be a B2FATAL, the module won't do anything useful if this happens
76 B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
77 std::vector <string> mod = {
78 "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
79 "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
80 "2.1.1", "2.1.2", "2.2.1", "2.2.2", "2.3.1", "2.3.2", "2.4.1", "2.4.2",
81 "2.5.1", "2.5.2", "2.6.1", "2.6.2", "2.7.1", "2.7.2", "2.8.1", "2.8.2",
82 "2.9.1", "2.9.2", "2.10.1", "2.10.2", "2.11.1", "2.11.2", "2.12.1", "2.12.2"
83 };
84 for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
85 // set some default size to nu, nv?
86 } else {
87 // Have been promised that all modules have the same number of pixels, so just take from the first one
89 nu = cellGetInfo.getUCells();
90 nv = cellGetInfo.getVCells();
91 }
92
93 for (VxdID& aPXDModule : m_PXDModules) {
94 auto buff = (std::string)aPXDModule;
95 replace(buff.begin(), buff.end(), '.', '_');
96 registerEpicsPV("PXD:Eff:" + buff, (std::string)aPXDModule);
97
98 TString histTitle = "PXD Hit Efficiency on Module " + (std::string)aPXDModule + ";Pixel in U;Pixel in V";
99 m_cEffModules[aPXDModule] = new TCanvas((m_histogramDirectoryName + "/c_Eff_" + buff).c_str());
100 m_eEffModules[aPXDModule] = new TEfficiency(("ePXDHitEff_" + buff).c_str(), histTitle,
101 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
102 }
103
104 m_cInnerMap = new TCanvas((m_histogramDirectoryName + "/c_InnerMap").data());
105 m_cOuterMap = new TCanvas((m_histogramDirectoryName + "/c_OuterMap").data());
106 m_hInnerMap = new TH2F("hEffInnerMap", "hEffInnerMap", m_u_bins * 8, 0, m_u_bins * 8, m_v_bins * 2, 0, m_v_bins * 2);
107 m_hOuterMap = new TH2F("hEffOuterMap", "hEffOuterMap", m_u_bins * 12, 0, m_u_bins * 12, m_v_bins * 2, 0, m_v_bins * 2);
108
109 m_nrxbins = m_PXDModules.size() + 3; // Modules + L1 + L2 + All
110 m_hErrorLine = new TH1F("hPXDErrorlimit", "Error Limit", m_nrxbins, 0, m_nrxbins);
111 m_hWarnLine = new TH1F("hPXDWarnlimit", "Warn Limit", m_nrxbins, 0, m_nrxbins);
112 for (int i = 0; i < (int)m_nrxbins; i++) {
113 m_hErrorLine->SetBinContent(i + 1, m_errorlevel);
114 m_hWarnLine->SetBinContent(i + 1, m_warnlevel);
115 }
116 m_hWarnLine->SetLineColor(kOrange - 3);
117 m_hWarnLine->SetLineWidth(3);
118 m_hWarnLine->SetLineStyle(4);
119 m_hErrorLine->SetLineColor(kRed + 3);
120 m_hErrorLine->SetLineWidth(3);
121 m_hErrorLine->SetLineStyle(7);
122
123 //One bin for each module in the geometry, one histogram for each layer
124 m_cEffAll = new TCanvas((m_histogramDirectoryName + "/c_EffAll").data());
125 m_eEffAll = new TEfficiency("ePXDHitEffAll", "PXD Integrated Efficiency of each module;PXD Module;", m_nrxbins, 0, m_nrxbins);
126 m_eEffAll->SetConfidenceLevel(m_confidence);
127 m_eEffAll->Paint("AP");
128 m_hEffAllLastTotal = m_eEffAll->GetCopyTotalHisto();
129 m_hEffAllLastPassed = m_eEffAll->GetCopyPassedHisto();
130
131 setLabels(m_eEffAll->GetPaintedGraph());
132
133 m_cEffAllUpdate = new TCanvas((m_histogramDirectoryName + "/c_EffAllUp").data());
134 m_eEffAllUpdate = new TEfficiency("ePXDHitEffAllUpdate", "PXD Integral and last-updated Efficiency per module;PXD Module;",
135 m_nrxbins, 0, m_nrxbins);
136 m_eEffAllUpdate->SetConfidenceLevel(m_confidence);
137
138 m_eEffAllUpdate->Paint("AP");
139 setLabels(m_eEffAllUpdate->GetPaintedGraph());
140
143
144 registerEpicsPV("PXD:Eff:Status", "Status");
145 registerEpicsPV("PXD:Eff:Overall", "All");
146 registerEpicsPV("PXD:Eff:L1", "L1");
147 registerEpicsPV("PXD:Eff:L2", "L2");
148 B2DEBUG(1, "DQMHistAnalysisPXDEff: initialized.");
149}
150
151
153{
154 B2DEBUG(1, "DQMHistAnalysisPXDEff: beginRun called.");
155
156 // Clear all used canvases
157 m_cEffAll->Clear();
158 m_cEffAllUpdate->Clear();
159 m_cInnerMap->Clear();
160 m_cOuterMap->Clear();
161 for (auto single_cmap : m_cEffModules) {
162 if (single_cmap.second) single_cmap.second->Clear();
163 }
164
165 // The 2d Efficiency maps (m_eEffModules[]) per module are not cleared, but re-created each update
166 // also they are only drawn to Canvas on update, thus no clear is needed here
167
168 // Reset TEfficiency and get (new) alarm limits from PVs
169 // no way to reset TEfficiency, do it bin by bin
170 for (int i = 0; i < m_nrxbins; i++) {
171 int bin = i + 1;
172 m_eEffAll->SetPassedEvents(bin, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
173 m_eEffAll->SetTotalEvents(bin, 0);
174 m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
175 m_eEffAllUpdate->SetTotalEvents(bin, 0);
176
177 if (i < int(m_PXDModules.size())) { // only for modules
180
181 // get warn and error limit
182 // as the same array as above, we assume chid exists
183 double dummy, loerr = 0, lowarn = 0;
184 if (requestLimitsFromEpicsPVs((std::string)m_PXDModules[i], loerr, lowarn, dummy, dummy)) {
185 m_hErrorLine->SetBinContent(bin, loerr);
187 m_hWarnLine->SetBinContent(bin, lowarn);
189 }
190 }
191 }
192 {
193 double dummy, loerr = 0, lowarn = 0;
196 if (requestLimitsFromEpicsPVs("L1", loerr, lowarn, dummy, dummy)) {
197 m_hErrorLine->SetBinContent(m_PXDModules.size() + 1, loerr);
198 if (m_perModuleAlarm) m_errorlevelmod["L1"] = loerr;
199 m_hWarnLine->SetBinContent(m_PXDModules.size() + 1, lowarn);
200 if (m_perModuleAlarm) m_warnlevelmod["L1"] = lowarn;
201 }
204 if (requestLimitsFromEpicsPVs("L2", loerr, lowarn, dummy, dummy)) {
205 m_hErrorLine->SetBinContent(m_PXDModules.size() + 2, loerr);
206 if (m_perModuleAlarm) m_errorlevelmod["L2"] = loerr;
207 m_hWarnLine->SetBinContent(m_PXDModules.size() + 2, lowarn);
208 if (m_perModuleAlarm) m_warnlevelmod["L2"] = lowarn;
209 }
212 if (requestLimitsFromEpicsPVs("All", loerr, lowarn, dummy, dummy)) {
213 m_hErrorLine->SetBinContent(m_PXDModules.size() + 3, loerr);
214 if (m_perModuleAlarm) m_errorlevelmod["All"] = loerr;
215 m_hWarnLine->SetBinContent(m_PXDModules.size() + 3, lowarn);
216 if (m_perModuleAlarm) m_warnlevelmod["All"] = lowarn;
217 }
218 }
219
220 // Clear all remaining Histograms (e.g. for our private delta histogramming)
221 m_hEffAllLastTotal->Reset();
222 m_hEffAllLastPassed->Reset();
223 m_hInnerMap->Reset();
224 m_hOuterMap->Reset();
225}
226
227bool DQMHistAnalysisPXDEffModule::updateEffBins(int bin, int nhit, int nmatch, int minentries)
228{
229 m_eEffAll->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
230 m_eEffAll->SetTotalEvents(bin, nhit);
231 m_eEffAll->SetPassedEvents(bin, nmatch);
232
233 if (nhit < minentries) {
234 // update the first entries directly (short runs)
235 m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
236 m_eEffAllUpdate->SetTotalEvents(bin, nhit);
237 m_eEffAllUpdate->SetPassedEvents(bin, nmatch);
238 m_hEffAllLastTotal->SetBinContent(bin, nhit);
239 m_hEffAllLastPassed->SetBinContent(bin, nmatch);
240 return true;
241 } else if (nhit - m_hEffAllLastTotal->GetBinContent(bin) > minentries) {
242 m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
243 m_eEffAllUpdate->SetTotalEvents(bin, nhit - m_hEffAllLastTotal->GetBinContent(bin));
244 m_eEffAllUpdate->SetPassedEvents(bin, nmatch - m_hEffAllLastPassed->GetBinContent(bin));
245 m_hEffAllLastTotal->SetBinContent(bin, nhit);
246 m_hEffAllLastPassed->SetBinContent(bin, nmatch);
247 return true;
248 }// else
249 return false;
250}
251
253{
254 bool warn_flag = (m_eEffAll->GetEfficiency(bin) + m_eEffAll->GetEfficiencyErrorUp(bin) <
255 m_warnlevelmod[name]); // (and not only the actual eff value)
256 if (m_alarmAdhoc) {
257 warn_flag |= (m_eEffAllUpdate->GetEfficiency(bin) + m_eEffAllUpdate->GetEfficiencyErrorUp(bin) <
258 m_warnlevelmod[name]); // (and not only the actual eff value)
259 }
260 return warn_flag;
261}
262
264{
265 bool error_flag = (m_eEffAll->GetEfficiency(bin) + m_eEffAll->GetEfficiencyErrorUp(bin) <
266 m_errorlevelmod[name]); // error if upper error value is below limit
267 if (m_alarmAdhoc) {
268 error_flag |= (m_eEffAllUpdate->GetEfficiency(bin) + m_eEffAllUpdate->GetEfficiencyErrorUp(bin) <
269 m_errorlevelmod[name]); // error if upper error value is below limit
270 }
271 return error_flag;
272}
273
274void DQMHistAnalysisPXDEffModule::setLabels(TGraphAsymmErrors* gr)
275{
276 if (gr) {
277 auto ax = gr->GetXaxis();
278 if (ax) {
279 ax->Set(m_nrxbins, 0, m_nrxbins);
280 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
281 TString ModuleName = (std::string)m_PXDModules[i];
282 ax->SetBinLabel(i + 1, ModuleName);
283 }
284 ax->SetBinLabel(m_PXDModules.size() + 1, "L1");
285 ax->SetBinLabel(m_PXDModules.size() + 2, "L2");
286 ax->SetBinLabel(m_PXDModules.size() + 3, "All");
287 }
288 }
289}
290
292{
293 {
294 // First create some 2d overview of efficiency for all modules
295 // This is not taken into account for efficiency calculation as
296 // there may be update glitches dues to separate histograms
297 // The histograms
298 bool updateinner = false, updateouter = false;
299 for (auto aPXDModule : m_PXDModules) {
300 auto buff = (std::string)aPXDModule;
301 replace(buff.begin(), buff.end(), '.', '_');
302
303 std::string locationHits = "track_hits_" + buff;
304 if (m_histogramDirectoryName != "") {
305 locationHits = m_histogramDirectoryName + "/" + locationHits;
306 }
307 std::string locationMatches = "matched_cluster_" + buff;
308 if (m_histogramDirectoryName != "") {
309 locationMatches = m_histogramDirectoryName + "/" + locationMatches;
310 }
311
312 auto Hits = findHist(locationHits, true);// check if updated
313 auto Matches = findHist(locationMatches, true);// check if updated
314
315 if (Hits == nullptr && Matches == nullptr) continue; // none updated
316
317 if (Hits == nullptr) Hits = findHist(locationHits); // actually, this should not happen ...
318 if (Matches == nullptr) Matches = findHist(locationMatches); // ... as updates should coincide
319
320 // Finding only one of them should only happen in very strange situations... still better check
321 if (Hits && Matches) {
322 if (m_cEffModules[aPXDModule] && m_eEffModules[aPXDModule]) {// this check creates them with a nullptr ..bad
323 m_eEffModules[aPXDModule]->SetTotalHistogram(*Hits, "f");
324 m_eEffModules[aPXDModule]->SetPassedHistogram(*Matches, "f");
325
326 m_cEffModules[aPXDModule]->cd();
327 m_eEffModules[aPXDModule]->Paint("colz"); // not Draw, enforce to create GetPaintedHistogram?
328 m_eEffModules[aPXDModule]->Draw("colz"); // but Draw needed to export Canvas!
329 m_cEffModules[aPXDModule]->Modified();
330 m_cEffModules[aPXDModule]->Update();
331 UpdateCanvas(m_cEffModules[aPXDModule]);
332
333 auto h = m_eEffModules[aPXDModule]->GetPaintedHistogram();
334 int s = (2 - aPXDModule.getSensorNumber()) * m_v_bins;
335 int l = (aPXDModule.getLadderNumber() - 1) * m_u_bins;
336 if (m_hInnerMap && aPXDModule.getLayerNumber() == 1) {
337 updateinner = true;
338 for (int u = 0; u < m_u_bins; u++) {
339 for (int v = 0; v < m_v_bins; v++) {
340 auto b = h->GetBin(u + 1, v + 1);
341 m_hInnerMap->Fill(u + l, v + s, h->GetBinContent(b));
342 }
343 }
344 }
345 if (m_hOuterMap && aPXDModule.getLayerNumber() == 2) {
346 updateouter = true;
347 for (int u = 0; u < m_u_bins; u++) {
348 for (int v = 0; v < m_v_bins; v++) {
349 auto b = h->GetBin(u + 1, v + 1);
350 m_hOuterMap->Fill(u + l, v + s, h->GetBinContent(b));
351 }
352 }
353 }
354 }
355 } else {
356 B2WARNING("only one plot upd " << aPXDModule);
357 }
358 }
359 // Single-Module histos + 2d overview finished. now draw overviews
360 if (updateinner) {
361 m_cInnerMap->cd();
362 if (m_hInnerMap) m_hInnerMap->Draw("colz");
363 m_cInnerMap->Modified();
364 m_cInnerMap->Update();
366 }
367 if (updateouter) {
368 m_cOuterMap->cd();
369 if (m_hOuterMap) m_hOuterMap->Draw("colz");
370 m_cOuterMap->Modified();
371 m_cOuterMap->Update();
373 }
374 // 3d overview done
375 }
376
377
378// Now, calculate and update efficiency.
379// Only histogram is used for numerator AND denominatorto to have an atomic update.
380// (avoid possible update glitches from daq/dqm framework side)
381// The bins per module are read out and filled into an TEfficiency as total and passed events into bin
382// Summaries for L1, L2, Overall are added, too
383 auto Combined = findHist(m_histogramDirectoryName, "PXD_Eff_combined", true);// only if updated
384
385 if (Combined) {
386 // only if histogram was changed
387
388 EStatus stat_data = c_StatusTooFew;
389 bool error_flag = false;
390 bool warn_flag = false;
391 double all = 0.0;
392
393 double imatch = 0.0, ihit = 0.0;
394 double imatchL1 = 0.0, ihitL1 = 0.0;
395 double imatchL2 = 0.0, ihitL2 = 0.0;
396 int ieff = 0; // count number of modules with useful stytistics
397
398 std::map <VxdID, bool> updated{}; // init to false, keep track of updated histograms
399 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
400 // workaround for excluded module
401 if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
402 // excluded modules are not counted at all!
403 int bin = i + 1; // bin nr is index +1
404
405 VxdID& aModule = m_PXDModules[i];
406 double nmatch = Combined->GetBinContent(i * 2 + 2);
407 double nhit = Combined->GetBinContent(i * 2 + 1);
408
409 imatch += nmatch;
410 ihit += nhit;
411 // check layer
412 if (i >= 16) {
413 imatchL2 += nmatch;
414 ihitL2 += nhit;
415 } else {
416 imatchL1 += nmatch;
417 ihitL1 += nhit;
418 }
419
420 if (nhit >= m_minEntries) { // dont update if there is nothing to calculate
421 ieff++; // only count in modules with significant stat
422 double var_e = nmatch / nhit; // can never be zero
423 m_monObj->setVariable(Form("efficiency_%d_%d_%d", aModule.getLayerNumber(), aModule.getLadderNumber(), aModule.getSensorNumber()),
424 var_e);
425 }
426
428 all += nhit;
429
430 updated[aModule] = updateEffBins(bin, nhit, nmatch, m_minEntries);
431
432 // workaround for excluded module
433 if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
434
435 // get the errors and check for limits for each bin separately ...
436
437 if (nhit >= m_minEntries) {
438 error_flag |= check_error_level(bin, aModule);
439 warn_flag |= check_warn_level(bin, aModule);
440 }
441 }
442
443 updateEffBins(m_PXDModules.size() + 1, ihitL1, imatchL1, m_minEntries * 8);
444 if (ihitL1 >= m_minEntries) {
445 error_flag |= check_error_level(m_PXDModules.size() + 1, "L1");
446 warn_flag |= check_warn_level(m_PXDModules.size() + 1, "L1");
447 }
448 updateEffBins(m_PXDModules.size() + 2, ihitL2, imatchL2, m_minEntries * 12);
449 if (ihitL2 >= m_minEntries) {
450 error_flag |= check_error_level(m_PXDModules.size() + 2, "L2");
451 warn_flag |= check_warn_level(m_PXDModules.size() + 2, "L2");
452 }
453 updateEffBins(m_PXDModules.size() + 3, ihit, imatch, m_minEntries * 20);
454 if (ihit >= m_minEntries) {
455 error_flag |= check_error_level(m_PXDModules.size() + 3, "All");
456 warn_flag |= check_warn_level(m_PXDModules.size() + 3, "All");
457 }
458
459 {
460 m_cEffAll->cd();
461 m_cEffAll->cd(0);
462 m_eEffAll->Paint("AP");
463 m_cEffAll->Clear();
464 m_cEffAll->cd(0);
465
466 auto gr = m_eEffAll->GetPaintedGraph();
467 if (gr) {
468 double scale_min = 1.0;
469 for (int i = 0; i < gr->GetN(); i++) {
470 gr->SetPointEXhigh(i, 0.);
471 gr->SetPointEXlow(i, 0.);
472 // this has to be done first, as it will recalc Min/Max and destroy axis
473 Double_t x, y;
474 gr->GetPoint(i, x, y);
475 gr->SetPoint(i, x - 0.01, y); // workaround for jsroot bug (fixed upstream)
476 auto val = y - gr->GetErrorYlow(i); // Error is relative to value
477 if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
478 // scale update only for included module
479 if (scale_min > val) scale_min = val;
480 }
481 }
482 if (scale_min == 1.0) scale_min = 0.0;
483 if (scale_min > 0.9) scale_min = 0.9;
484 auto ay = gr->GetYaxis();
485 if (ay) ay->SetRangeUser(scale_min, 1.0);
486 setLabels(gr);
487
488 gr->SetLineColor(4);
489 gr->SetLineWidth(2);
490 gr->SetMarkerStyle(8);
491
492 gr->Draw("AP");
493
494 for (auto& it : m_excluded) {
495 static std::map <int, TLatex*> ltmap;
496 auto tt = ltmap[it];
497 if (!tt) {
498 tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
499 tt->SetTextSize(0.035);
500 tt->SetTextAngle(90);// Rotated
501 tt->SetTextAlign(12);// Centered
502 ltmap[it] = tt;
503 } else {
504 tt->SetY(scale_min);
505 }
506 tt->Draw();
507 }
508
509
510 EStatus all_stat = makeStatus(all >= m_minEntries, warn_flag, error_flag);
511 colorizeCanvas(m_cEffAll, all_stat);
512
513 m_hWarnLine->Draw("same,hist");
514 m_hErrorLine->Draw("same,hist");
515 }
516
517 UpdateCanvas(m_cEffAll->GetName());
518 m_cEffAll->Modified();
519 m_cEffAll->Update();
520 }
521
522 {
523 m_cEffAllUpdate->cd();
524 m_eEffAllUpdate->Paint("AP");
525 m_cEffAllUpdate->Clear();
526 m_cEffAllUpdate->cd(0);
527
528 auto gr = m_eEffAllUpdate->GetPaintedGraph();
529 auto gr3 = (TGraphAsymmErrors*) m_eEffAll->GetPaintedGraph()->Clone();
530 if (gr3) {
531 for (int i = 0; i < gr3->GetN(); i++) {
532 Double_t x, y;
533 gr3->GetPoint(i, x, y);
534 gr3->SetPoint(i, x + 0.2, y);
535 }
536 }
537
538 double scale_min = 1.0;
539 if (gr) {
540 for (int i = 0; i < gr->GetN(); i++) {
541 gr->SetPointEXhigh(i, 0.);
542 gr->SetPointEXlow(i, 0.);
543 // this has to be done first, as it will recalc Min/Max and destroy axis
544 Double_t x, y;
545 gr->GetPoint(i, x, y);
546 gr->SetPoint(i, x - 0.2, y); // shift a bit if in same plot
547 auto val = y - gr->GetErrorYlow(i); // Error is relative to value
548 if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
549 // skip scale update only for included modules
550 if (scale_min > val) scale_min = val;
551 }
552 }
553 if (scale_min == 1.0) scale_min = 0.0;
554 if (scale_min > 0.9) scale_min = 0.9;
555 auto ay = gr->GetYaxis();
556 if (ay) ay->SetRangeUser(scale_min, 1.0);
557 setLabels(gr);
558
559 for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
560 if (updated[m_PXDModules[i]]) {
561 // we should only write if it was updated!
562 Double_t x, y;// we assume that double and Double_t are same!
563 gr->GetPoint(i, x, y);
564 setEpicsPV((std::string)m_PXDModules[i], y);
565 }
566 }
567
568 gr->SetLineColor(kBlack);
569 gr->SetLineWidth(3);
570 gr->SetMarkerStyle(33);
571 gr->Draw("AP");
572 } else scale_min = 0.0;
573 if (gr3) gr3->Draw("P"); // both in one plot
574
575 for (auto& it : m_excluded) {
576 std::map <int, TLatex*> ltmap;
577 auto tt = ltmap[it];
578 if (!tt) {
579 tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
580 tt->SetTextSize(0.035);
581 tt->SetTextAngle(90);// Rotated
582 tt->SetTextAlign(12);// Centered
583 } else {
584 tt->SetY(scale_min);
585 }
586 tt->Draw();
587 }
588
589 stat_data = makeStatus(all >= m_minEntries, warn_flag, error_flag);
591
592 m_hWarnLine->Draw("same,hist");
593 m_hErrorLine->Draw("same,hist");
594 }
595 UpdateCanvas(m_cEffAllUpdate->GetName());
596 m_cEffAllUpdate->Modified();
597 m_cEffAllUpdate->Update();
598
599 double var_efficiency = ihit > 0 ? imatch / ihit : 0.0;
600 double var_efficiencyL1 = ihitL1 > 0 ? imatchL1 / ihitL1 : 0.0;
601 double var_efficiencyL2 = ihitL2 > 0 ? imatchL2 / ihitL2 : 0.0;
602
603 m_monObj->setVariable("efficiency", var_efficiency);
604 m_monObj->setVariable("efficiencyL1", var_efficiencyL1);
605 m_monObj->setVariable("efficiencyL2", var_efficiencyL2);
606 m_monObj->setVariable("nmodules", ieff);
607
608 setEpicsPV("Status", stat_data);
609 // only update if statistics is reasonable, we dont want "0" drops between runs!
610 if (stat_data != c_StatusTooFew) {
611 setEpicsPV("All", var_efficiency);
612 setEpicsPV("L1", var_efficiencyL1);
613 setEpicsPV("L2", var_efficiencyL2);
614 }
615 }
616}
617
619{
620 B2DEBUG(1, "DQMHistAnalysisPXDEff: terminate called");
621
622 for (VxdID& aPXDModule : m_PXDModules) {
623 if (m_cEffModules[aPXDModule]) delete m_cEffModules[aPXDModule];
624 if (m_eEffModules[aPXDModule]) delete m_eEffModules[aPXDModule];
625 }
626
629
630 if (m_cInnerMap) delete m_cInnerMap;
631 if (m_cOuterMap) delete m_cOuterMap;
632 if (m_hInnerMap) delete m_hInnerMap;
633 if (m_hOuterMap) delete m_hOuterMap;
634
635 if (m_hErrorLine) delete m_hErrorLine;
636 if (m_hWarnLine) delete m_hWarnLine;
637
638 if (m_cEffAll) delete m_cEffAll;
639 if (m_eEffAll) delete m_eEffAll;
640
643}
644
The base class for the histogram analysis module.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
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).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
EStatus
Status flag of histogram/canvas.
@ c_StatusTooFew
Not enough entries/event to judge.
EStatus makeStatus(bool enough, bool warn_flag, bool error_flag)
Helper function to judge the status for coloring and EPICS.
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.
TCanvas * m_cOuterMap
Full Eff Map Outer Layer.
TCanvas * m_cInnerMap
Full Eff Map Inner Layer.
void terminate(void) override final
This method is called at the end of the event processing.
std::map< VxdID, TEfficiency * > m_eEffModules
Individual efficiency for each module, 2d histogram.
bool m_perModuleAlarm
use alarm level per module
std::map< VxdID, TCanvas * > m_cEffModules
Individual efficiency for each module, canvas.
void initialize(void) override final
Initializer.
int m_nrxbins
Number of bins in efficiency plot, all modules plus layer and summary.
TH1 * m_hEffAllLastTotal
TH1, last state, total.
TH1 * m_hEffAllLastPassed
TH1, last state, passed.
TH2F * m_hOuterMap
Full Eff Map Outer Layer.
void setLabels(TGraphAsymmErrors *gr)
Set module labels for TGraphAsymmErrors.
TEfficiency * m_eEffAll
One bin for each module in the geometry.
MonitoringObject * m_monObj
Monitoring Object.
bool check_warn_level(int bin, std::string name)
Check bin/name for warn condition.
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
std::string m_histogramDirectoryName
name of histogram directory
TCanvas * m_cEffAllUpdate
Final Canvas for Update.
TH2F * m_hInnerMap
Full Eff Map Inner Layer.
double m_confidence
confidence level for error bars
std::map< std::string, double > m_warnlevelmod
warn level for alarm per module
TEfficiency * m_eEffAllUpdate
Efficiency, last state, updated.
TH1F * m_hWarnLine
TLine object for warning limit.
double m_errorlevel
error level for alarm
bool check_error_level(int bin, std::string name)
Check bin/name for error condition.
std::map< std::string, double > m_errorlevelmod
error level for alarm per module
TH1F * m_hErrorLine
TLine object for error error.
std::vector< int > m_excluded
Indizes of excluded PXD Modules.
bool m_alarmAdhoc
generate alarm from adhoc values
bool updateEffBins(int bin, int nhit, int nmatch, int minentries)
Update bin in efficiency plots with condition on nhits.
double m_warnlevel
warn level for alarm
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
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 addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
int getVCells() const
Return number of pixel/strips in v direction.
int getUCells() const
Return number of pixel/strips in u direction.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
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.
STL namespace.