Belle II Software  release-08-00-10
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 
14 #include <dqm/analysis/modules/DQMHistAnalysisPXDEff.h>
15 #include <TROOT.h>
16 #include <TLatex.h>
17 #include <TGraphAsymmErrors.h>
18 #include <vxd/geometry/GeoCache.h>
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 //-----------------------------------------------------------------
24 // Register the Module
25 //-----------------------------------------------------------------
26 REG_MODULE(DQMHistAnalysisPXDEff);
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
32 DQMHistAnalysisPXDEffModule::DQMHistAnalysisPXDEffModule() : DQMHistAnalysisModule()
33 {
34  // This module CAN NOT be run in parallel!
35  setDescription("DQM Analysis for PXD Efficiency");
36 
37  // Parameter definition
38 
39  // 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
40  // Or get one and clone it
41  addParam("binsU", m_u_bins, "histogram bins in u direction, needs to be the same as in PXDDQMEfficiency", int(16));
42  addParam("binsV", m_v_bins, "histogram bins in v direction, needs to be the same as in PXDDQMEfficiency", int(48));
43  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms were placed",
44  std::string("PXDEFF"));
45  addParam("ConfidenceLevel", m_confidence, "Confidence Level for error bars and alarms", 0.9544);
46  addParam("WarnLevel", m_warnlevel, "Efficiency Warn Level for alarms", 0.92);
47  addParam("ErrorLevel", m_errorlevel, "Efficiency Level for alarms", 0.90);
48  addParam("perModuleAlarm", m_perModuleAlarm, "Alarm level per module", true);
49  addParam("alarmAdhoc", m_alarmAdhoc, "Generate Alarm from adhoc values", true);
50  addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 1000);
51  addParam("excluded", m_excluded, "excluded module (indizes starting from 0 to 39)");
52  B2DEBUG(1, "DQMHistAnalysisPXDEff: Constructor done.");
53 }
54 
56 {
57 }
58 
60 {
61  B2DEBUG(99, "DQMHistAnalysisPXDEffModule: initialized.");
62 
65 
66  // collect the list of all PXD Modules in the geometry here
67  std::vector<VxdID> sensors = geo.getListOfSensors();
68  for (VxdID& aVxdID : sensors) {
69  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
70  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
71  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
72  }
73  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
74 
75  gROOT->cd(); // this seems to be important, or strange things happen
76 
77  int nu = 1;//If this does not get overwritten, the histograms will anyway never contain anything useful
78  int nv = 1;
79  if (m_PXDModules.size() == 0) {
80  // This could as well be a B2FATAL, the module won't do anything useful if this happens
81  B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
82  std::vector <string> mod = {
83  "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
84  "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
85  "2.1.1", "2.1.2", "2.2.1", "2.2.2", "2.3.1", "2.3.2", "2.4.1", "2.4.2",
86  "2.5.1", "2.5.2", "2.6.1", "2.6.2", "2.7.1", "2.7.2", "2.8.1", "2.8.2",
87  "2.9.1", "2.9.2", "2.10.1", "2.10.2", "2.11.1", "2.11.2", "2.12.1", "2.12.2"
88  };
89  for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
90  // set some default size to nu, nv?
91  } else {
92  // Have been promised that all modules have the same number of pixels, so just take from the first one
93  VXD::SensorInfoBase cellGetInfo = geo.getSensorInfo(m_PXDModules[0]);
94  nu = cellGetInfo.getUCells();
95  nv = cellGetInfo.getVCells();
96  }
97 
98  for (VxdID& aPXDModule : m_PXDModules) {
99  auto buff = (std::string)aPXDModule;
100  replace(buff.begin(), buff.end(), '.', '_');
101  registerEpicsPV("PXD:Eff:" + buff, (std::string)aPXDModule);
102 
103  TString histTitle = "PXD Hit Efficiency on Module " + (std::string)aPXDModule + ";Pixel in U;Pixel in V";
104  m_cEffModules[aPXDModule] = new TCanvas((m_histogramDirectoryName + "/c_Eff_" + buff).c_str());
105  m_hEffModules[aPXDModule] = new TEfficiency(("ePXDHitEff_" + buff).c_str(), histTitle,
106  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
107  }
108 
109  m_cInnerMap = new TCanvas((m_histogramDirectoryName + "/c_InnerMap").data());
110  m_cOuterMap = new TCanvas((m_histogramDirectoryName + "/c_OuterMap").data());
111  m_hInnerMap = new TH2F("hEffInnerMap", "hEffInnerMap", m_u_bins * 8, 0, m_u_bins * 8, m_v_bins * 2, 0, m_v_bins * 2);
112  m_hOuterMap = new TH2F("hEffOuterMap", "hEffOuterMap", m_u_bins * 12, 0, m_u_bins * 12, m_v_bins * 2, 0, m_v_bins * 2);
113 
114  m_hErrorLine = new TH1F("hPXDErrorlimit", "Error Limit", m_PXDModules.size(), 0, m_PXDModules.size());
115  m_hWarnLine = new TH1F("hPXDWarnlimit", "Warn Limit", m_PXDModules.size(), 0, m_PXDModules.size());
116  for (int i = 0; i < (int)m_PXDModules.size(); i++) {
117  m_hErrorLine->SetBinContent(i + 1, m_errorlevel);
118  m_hWarnLine->SetBinContent(i + 1, m_warnlevel);
119  }
120  m_hWarnLine->SetLineColor(kOrange - 3);
121  m_hWarnLine->SetLineWidth(3);
122  m_hWarnLine->SetLineStyle(4);
123  m_hErrorLine->SetLineColor(kRed + 3);
124  m_hErrorLine->SetLineWidth(3);
125  m_hErrorLine->SetLineStyle(7);
126 
127  //One bin for each module in the geometry, one histogram for each layer
128  m_cEffAll = new TCanvas((m_histogramDirectoryName + "/c_EffAll").data());
129  m_hEffAll = new TEfficiency("ePXDHitEffAll", "PXD Integrated Efficiency of each module;PXD Module;",
130  m_PXDModules.size(), 0, m_PXDModules.size());
131  m_hEffAll->SetConfidenceLevel(m_confidence);
132  m_hEffAll->Paint("AP");
133  m_hEffAllLastTotal = m_hEffAll->GetCopyTotalHisto();
134  m_hEffAllLastPassed = m_hEffAll->GetCopyPassedHisto();
135 
136  {
137  auto gr = m_hEffAll->GetPaintedGraph();
138 
139  if (gr) {
140  auto ax = gr->GetXaxis();
141  if (ax) {
142  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
143  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
144  TString ModuleName = (std::string)m_PXDModules[i];
145  ax->SetBinLabel(i + 1, ModuleName);
146  }
147  }
148  }
149  }
150 
151  m_cEffAllUpdate = new TCanvas((m_histogramDirectoryName + "/c_EffAllUp").data());
152  m_hEffAllUpdate = new TEfficiency("ePXDHitEffAllUpdate", "PXD Integral and last-updated Efficiency per module;PXD Module;",
153  m_PXDModules.size(), 0, m_PXDModules.size());
154  m_hEffAllUpdate->SetConfidenceLevel(m_confidence);
155 
156  {
157  auto gr = m_hEffAllUpdate->GetPaintedGraph();
158 
159  if (gr) {
160  auto ax = gr->GetXaxis();
161  if (ax) {
162  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
163  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
164  TString ModuleName = (std::string)m_PXDModules[i];
165  ax->SetBinLabel(i + 1, ModuleName);
166  }
167  }
168  }
169  }
170 
173 
174  registerEpicsPV("PXD:Eff:Status", "Status");
175  registerEpicsPV("PXD:Eff:Overall", "Overall");
176  B2DEBUG(1, "DQMHistAnalysisPXDEff: initialized.");
177 }
178 
179 
181 {
182  B2DEBUG(1, "DQMHistAnalysisPXDEff: beginRun called.");
183 
184  m_cEffAll->Clear();
185  m_cEffAllUpdate->Clear();
186 
187  // no way to reset TEfficiency, do it bin by bin
188  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
189  int j = i + 1;
190  m_hEffAll->SetPassedEvents(j, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
191  m_hEffAll->SetTotalEvents(j, 0);
192  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
193  m_hEffAllUpdate->SetTotalEvents(j, 0);
194 
197 
198  // get warn and error limit
199  // as the same array as above, we assume chid exists
200  double dummy, loerr = 0, lowarn = 0;
201  if (requestLimitsFromEpicsPVs((std::string)m_PXDModules[i], loerr, lowarn, dummy, dummy)) {
202  m_hErrorLine->SetBinContent(i + 1, loerr);
204  m_hWarnLine->SetBinContent(i + 1, lowarn);
205  if (m_perModuleAlarm) m_warnlevelmod[m_PXDModules[i]] = lowarn;
206  }
207 
208  }
209  // Thus histo will contain old content until first update
210  m_hEffAllLastTotal->Reset();
211  m_hEffAllLastPassed->Reset();
212 
213  for (auto single_cmap : m_cEffModules) {
214  if (single_cmap.second) single_cmap.second->Clear();
215  }
216 
217 }
218 
219 
221 {
222  //Count how many of each type of histogram there are for the averaging
223  //std::map<std::string, int> typeCounter;
224 
225  m_hInnerMap->Reset();
226  m_hOuterMap->Reset();
227 
228  for (unsigned int i = 1; i <= m_PXDModules.size(); i++) {
229  VxdID& aPXDModule = m_PXDModules[i - 1];
230 
231  TString buff = (std::string)aPXDModule;
232  buff.ReplaceAll(".", "_");
233 
234  TH1* Hits, *Matches;
235  TString locationHits = "track_hits_" + buff;
236  if (m_histogramDirectoryName != "") {
237  locationHits = m_histogramDirectoryName + "/" + locationHits;
238  }
239  Hits = (TH1*)findHist(locationHits.Data());
240  TString locationMatches = "matched_cluster_" + buff;
241  if (m_histogramDirectoryName != "") {
242  locationMatches = m_histogramDirectoryName + "/" + locationMatches;
243  }
244  Matches = (TH1*)findHist(locationMatches.Data());
245 
246  // Finding only one of them should only happen in very strange situations...
247  if (Hits == nullptr || Matches == nullptr) {
248  B2ERROR("Missing histogram for sensor " << aPXDModule);
249  } else {
250  if (m_cEffModules[aPXDModule] && m_hEffModules[aPXDModule]) {// this check creates them with a nullptr ..bad
251  m_hEffModules[aPXDModule]->SetTotalHistogram(*Hits, "f");
252  m_hEffModules[aPXDModule]->SetPassedHistogram(*Matches, "f");
253 
254  m_cEffModules[aPXDModule]->cd();
255  m_hEffModules[aPXDModule]->Paint("colz"); // not Draw, enforce to create GetPaintedHistogram?
256  m_cEffModules[aPXDModule]->Modified();
257  m_cEffModules[aPXDModule]->Update();
258 
259  auto h = m_hEffModules[aPXDModule]->GetPaintedHistogram();
260  int s = (2 - aPXDModule.getSensorNumber()) * m_v_bins;
261  int l = (aPXDModule.getLadderNumber() - 1) * m_u_bins;
262  if (m_hInnerMap && aPXDModule.getLayerNumber() == 1) {
263  for (int u = 0; u < m_u_bins; u++) {
264  for (int v = 0; v < m_v_bins; v++) {
265  auto b = h->GetBin(u + 1, v + 1);
266  m_hInnerMap->Fill(u + l, v + s, h->GetBinContent(b));
267  }
268  }
269  }
270  if (m_hOuterMap && aPXDModule.getLayerNumber() == 2) {
271  for (int u = 0; u < m_u_bins; u++) {
272  for (int v = 0; v < m_v_bins; v++) {
273  auto b = h->GetBin(u + 1, v + 1);
274  m_hOuterMap->Fill(u + l, v + s, h->GetBinContent(b));
275  }
276  }
277  }
278  }
279  }
280  }// Single-Module histos + 2d overview finished
281 
282  m_cInnerMap->cd();
283  if (m_hInnerMap) m_hInnerMap->Draw("colz");
284  m_cInnerMap->Modified();
285  m_cInnerMap->Update();
286  m_cOuterMap->cd();
287  if (m_hOuterMap) m_hOuterMap->Draw("colz");
288  m_cOuterMap->Modified();
289  m_cOuterMap->Update();
290  UpdateCanvas(m_cInnerMap->GetName());
291  UpdateCanvas(m_cOuterMap->GetName());
292 
293  // Change: We now use one histogram for hits and matches to make
294  // sure that we have an atomic update which is otherwise not
295  // guaranteed by DQM framework
296  TString locationHits = "PXD_Eff_combined";
297  if (m_histogramDirectoryName != "") {
298  locationHits = m_histogramDirectoryName + "/" + locationHits;
299  }
300  TH1* Combined = (TH1*)findHist(locationHits.Data());
301 
302  double stat_data = 0;
303  bool error_flag = false;
304  bool warn_flag = false;
305  double all = 0.0;
306 
307  double imatch = 0.0, ihit = 0.0;
308  int ieff = 0;
309 
310  std::map <VxdID, bool> updated{}; // init to false, keep track of updated histograms
311  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
312  int j = i + 1;
313 
314  if (Combined == nullptr) {
315  m_hEffAll->SetPassedEvents(j, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
316  m_hEffAll->SetTotalEvents(j, 0);
317  } else {
318  VxdID& aModule = m_PXDModules[i];
319  double nmatch = Combined->GetBinContent(i * 2 + 2);
320  double nhit = Combined->GetBinContent(i * 2 + 1);
321  if (nmatch > 10 && nhit > 10) { // could be zero, too
322  imatch += nmatch;
323  ihit += nhit;
324  ieff++; // only count in modules working
325  double var_e = nmatch / nhit; // can never be zero
326  // workaround for excluded module
327  if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
328 
329  m_monObj->setVariable(Form("efficiency_%d_%d_%d", aModule.getLayerNumber(), aModule.getLadderNumber(), aModule.getSensorNumber()),
330  var_e);
331  }
332 
334  all += nhit;
335  m_hEffAll->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
336  m_hEffAll->SetTotalEvents(j, nhit);
337  m_hEffAll->SetPassedEvents(j, nmatch);
338 
339  if (nhit < m_minEntries) {
340  // update the first entries directly (short runs)
341  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
342  m_hEffAllUpdate->SetTotalEvents(j, nhit);
343  m_hEffAllUpdate->SetPassedEvents(j, nmatch);
344  m_hEffAllLastTotal->SetBinContent(j, nhit);
345  m_hEffAllLastPassed->SetBinContent(j, nmatch);
346  updated[aModule] = true;
347  } else if (nhit - m_hEffAllLastTotal->GetBinContent(j) > m_minEntries) {
348  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
349  m_hEffAllUpdate->SetTotalEvents(j, nhit - m_hEffAllLastTotal->GetBinContent(j));
350  m_hEffAllUpdate->SetPassedEvents(j, nmatch - m_hEffAllLastPassed->GetBinContent(j));
351  m_hEffAllLastTotal->SetBinContent(j, nhit);
352  m_hEffAllLastPassed->SetBinContent(j, nmatch);
353  updated[aModule] = true;
354  }
355 
356  // workaround for excluded module
357  if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
358 
359  // get the errors and check for limits for each bin seperately ...
360 
361  if (nhit > 50) {
362  error_flag |= (m_hEffAll->GetEfficiency(j) + m_hEffAll->GetEfficiencyErrorUp(j) <
363  m_errorlevelmod[aModule]); // error if upper error value is below limit
364  warn_flag |= (m_hEffAll->GetEfficiency(j) + m_hEffAll->GetEfficiencyErrorUp(j) <
365  m_warnlevelmod[aModule]); // (and not only the actual eff value)
366  if (m_alarmAdhoc) {
367  error_flag |= (m_hEffAllUpdate->GetEfficiency(j) + m_hEffAllUpdate->GetEfficiencyErrorUp(j) <
368  m_errorlevelmod[aModule]); // error if upper error value is below limit
369  warn_flag |= (m_hEffAllUpdate->GetEfficiency(j) + m_hEffAllUpdate->GetEfficiencyErrorUp(j) <
370  m_warnlevelmod[aModule]); // (and not only the actual eff value)
371  }
372  }
373  }
374  }
375 
376  {
377  m_cEffAll->cd();
378  m_cEffAll->cd(0);
379  m_hEffAll->Paint("AP");
380  m_cEffAll->Clear();
381  m_cEffAll->cd(0);
382 
383  auto gr = m_hEffAll->GetPaintedGraph();
384  if (gr) {
385  double scale_min = 1.0;
386  for (int i = 0; i < gr->GetN(); i++) {
387  gr->SetPointEXhigh(i, 0.);
388  gr->SetPointEXlow(i, 0.);
389  // this has to be done first, as it will recalc Min/Max and destroy axis
390  Double_t x, y;
391  gr->GetPoint(i, x, y);
392  gr->SetPoint(i, x - 0.01, y); // workaround for jsroot bug (fixed upstream)
393  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
394  if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
395  // scale update only for included module
396  if (scale_min > val) scale_min = val;
397  }
398  }
399  if (scale_min == 1.0) scale_min = 0.0;
400  if (scale_min > 0.9) scale_min = 0.9;
401  gr->SetMinimum(0);
402  gr->SetMaximum(m_PXDModules.size());
403  auto ay = gr->GetYaxis();
404  if (ay) ay->SetRangeUser(scale_min, 1.0);
405  auto ax = gr->GetXaxis();
406  if (ax) {
407  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
408  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
409  TString ModuleName = (std::string)m_PXDModules[i];
410  ax->SetBinLabel(i + 1, ModuleName);
411  }
412  }
413 
414  gr->SetLineColor(4);
415  gr->SetLineWidth(2);
416  gr->SetMarkerStyle(8);
417 
418  gr->Draw("AP");
419 
420  for (auto& it : m_excluded) {
421  auto tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
422  tt->SetTextAngle(90);// Rotated
423  tt->SetTextAlign(12);// Centered
424  tt->Draw();
425  }
426 
427  if (all < 100.) {
428  m_cEffAll->Pad()->SetFillColor(kGray);// Magenta or Gray
429  } else {
430  if (error_flag) {
431  m_cEffAll->Pad()->SetFillColor(kRed);// Red
432  } else if (warn_flag) {
433  m_cEffAll->Pad()->SetFillColor(kYellow);// Yellow
434  } else {
435  m_cEffAll->Pad()->SetFillColor(kGreen);// Green
436  // m_cEffAll->Pad()->SetFillColor(kWhite);// White
437  }
438  }
439 
440  m_cEffAll->Pad()->SetFrameFillColor(kWhite - 1); // White
441  m_cEffAll->Pad()->SetFrameFillStyle(1001);// White
442  m_cEffAll->Pad()->Modified();
443  m_cEffAll->Pad()->Update();
444  m_hWarnLine->Draw("same,hist");
445  m_hErrorLine->Draw("same,hist");
446  }
447 
448  UpdateCanvas(m_cEffAll->GetName());
449  m_cEffAll->Modified();
450  m_cEffAll->Update();
451  }
452 
453  {
454  m_cEffAllUpdate->cd();
455  m_hEffAllUpdate->Paint("AP");
456  m_cEffAllUpdate->Clear();
457  m_cEffAllUpdate->cd(0);
458 
459  auto gr = m_hEffAllUpdate->GetPaintedGraph();
460  auto gr3 = (TGraphAsymmErrors*) m_hEffAll->GetPaintedGraph()->Clone();
461  if (gr3) {
462  for (int i = 0; i < gr3->GetN(); i++) {
463  Double_t x, y;
464  gr3->GetPoint(i, x, y);
465  gr3->SetPoint(i, x + 0.2, y);
466  }
467  }
468 
469  double scale_min = 1.0;
470  if (gr) {
471  for (int i = 0; i < gr->GetN(); i++) {
472  gr->SetPointEXhigh(i, 0.);
473  gr->SetPointEXlow(i, 0.);
474  // this has to be done first, as it will recalc Min/Max and destroy axis
475  Double_t x, y;
476  gr->GetPoint(i, x, y);
477  gr->SetPoint(i, x - 0.2, y); // shift a bit if in same plot
478  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
479  if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
480  // skip scale update only for included modules
481  if (scale_min > val) scale_min = val;
482  }
483  }
484  if (scale_min == 1.0) scale_min = 0.0;
485  if (scale_min > 0.9) scale_min = 0.9;
486  gr->SetMinimum(0);
487  gr->SetMaximum(m_PXDModules.size());
488  auto ay = gr->GetYaxis();
489  if (ay) ay->SetRangeUser(scale_min, 1.0);
490  auto ax = gr->GetXaxis();
491  if (ax) {
492  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
493  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
494  TString ModuleName = (std::string)m_PXDModules[i];
495  ax->SetBinLabel(i + 1, ModuleName);
496  }
497  }
498  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
499  if (updated[m_PXDModules[i]]) {
500  // we should only write if it was updated!
501  Double_t x, y;// we assume that double and Double_t are same!
502  gr->GetPoint(i, x, y);
503  setEpicsPV((std::string)m_PXDModules[i], y);
504  }
505  }
506 
507  gr->SetLineColor(kBlack);
508  gr->SetLineWidth(3);
509  gr->SetMarkerStyle(33);
510  gr->Draw("AP");
511  } else scale_min = 0.0;
512  if (gr3) gr3->Draw("P"); // both in one plot
513 
514  for (auto& it : m_excluded) {
515  auto tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
516  tt->SetTextSize(0.035);
517  tt->SetTextAngle(90);// Rotated
518  tt->SetTextAlign(12);// Centered
519  tt->Draw();
520  }
521 
522  if (all < 100.) {
523  m_cEffAllUpdate->Pad()->SetFillColor(kGray);// Magenta or Gray
524  stat_data = 0.;
525  } else {
526  if (error_flag) {
527  m_cEffAllUpdate->Pad()->SetFillColor(kRed);// Red
528  stat_data = 4.;
529  } else if (warn_flag) {
530  m_cEffAllUpdate->Pad()->SetFillColor(kYellow);// Yellow
531  stat_data = 3.;
532  } else {
533  m_cEffAllUpdate->Pad()->SetFillColor(kGreen);// Green
534  stat_data = 2.;
536  // m_cEffAllUpdate->Pad()->SetFillColor(kWhite);// White
537  }
538  }
539  m_cEffAllUpdate->Pad()->SetFrameFillColor(kWhite - 1); // White
540  m_cEffAllUpdate->Pad()->SetFrameFillStyle(1001);// White
541  m_cEffAllUpdate->Pad()->Modified();
542  m_cEffAllUpdate->Pad()->Update();
543  m_hWarnLine->Draw("same,hist");
544  m_hErrorLine->Draw("same,hist");
545  }
546  UpdateCanvas(m_cEffAllUpdate->GetName());
547  m_cEffAllUpdate->Modified();
548  m_cEffAllUpdate->Update();
549 
550 
551  double var_efficiency = ihit > 0 ? imatch / ihit : 0.0;
552  m_monObj->setVariable("efficiency", var_efficiency);
553  m_monObj->setVariable("nmodules", ieff);
554 
555  setEpicsPV("Status", stat_data);
556  // only update if statistics is reasonable, we dont want "0" drops between runs!
557  if (stat_data != 0) setEpicsPV("Overall", var_efficiency);
558 }
559 
561 {
562  B2DEBUG(1, "DQMHistAnalysisPXDEff: terminate called");
563 }
564 
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
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.
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.
int m_minEntries
Update entry intervall.
std::map< VxdID, double > m_warnlevelmod
warn level for alarm per module
bool m_perModuleAlarm
use alarm level per module
TEfficiency * m_hEffAll
One bin for each module in the geometry.
std::map< VxdID, TCanvas * > m_cEffModules
Individual efficiency for each module, canvas.
void initialize(void) override final
Initializer.
TEfficiency * m_hEffAllUpdate
Efficiency, last state, updated.
std::map< VxdID, double > m_errorlevelmod
error level for alarm per module
TH1 * m_hEffAllLastTotal
TH1, last state, total.
TH1 * m_hEffAllLastPassed
TH1, last state, passed.
TH2F * m_hOuterMap
Full Eff Map Outer Layer.
MonitoringObject * m_monObj
Monitoring Object.
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
TH1F * m_hWarnLine
TLine object for warning limit.
double m_errorlevel
error level for alarm
std::map< VxdID, TEfficiency * > m_hEffModules
Individual efficiency for each module, 2d histogram.
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
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.