Belle II Software  release-06-02-00
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 
33 {
34  // This module CAN NOT be run in parallel!
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("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:Eff:"));
45  addParam("useEpics", m_useEpics, "Whether to update EPICS PVs.", false);
46  addParam("useEpicsRO", m_useEpicsRO, "useEpics ReadOnly", false);
47  addParam("ConfidenceLevel", m_confidence, "Confidence Level for error bars and alarms", 0.9544);
48  addParam("WarnLevel", m_warnlevel, "Efficiency Warn Level for alarms", 0.92);
49  addParam("ErrorLevel", m_errorlevel, "Efficiency Level for alarms", 0.90);
50  addParam("perModuleAlarm", m_perModuleAlarm, "Alarm level per module", true);
51  addParam("alarmAdhoc", m_alarmAdhoc, "Generate Alarm from adhoc values", true);
52  addParam("minEntries", m_minEntries, "minimum number of new entries for last time slot", 1000);
53  B2DEBUG(1, "DQMHistAnalysisPXDEff: Constructor done.");
54 }
55 
56 DQMHistAnalysisPXDEffModule::~DQMHistAnalysisPXDEffModule()
57 {
58 #ifdef _BELLE2_EPICS
59  if (m_useEpics) {
60  if (ca_current_context()) ca_context_destroy();
61  }
62 #endif
63 }
64 
65 void DQMHistAnalysisPXDEffModule::initialize()
66 {
67  B2DEBUG(99, "DQMHistAnalysisPXDEffModule: initialized.");
68 
69  m_monObj = getMonitoringObject("pxd");
70  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
71 
72  // collect the list of all PXD Modules in the geometry here
73  std::vector<VxdID> sensors = geo.getListOfSensors();
74  for (VxdID& aVxdID : sensors) {
75  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
76  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
77  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
78  }
79  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
80 
81  gROOT->cd(); // this seems to be important, or strange things happen
82 
83  int nu = 1;//If this does not get overwritten, the histograms will anyway never contain anything useful
84  int nv = 1;
85  //Have been promised that all modules have the same number of pixels, so just take from the first one
86  if (m_PXDModules.size() == 0) {
87  //This could as well be a B2FATAL, the module won't do anything useful if this happens
88  B2WARNING("No PXDModules in Geometry found! Use hard-coded setup.");
89  std::vector <string> mod = {
90  "1.1.1", "1.1.2", "1.2.1", "1.2.2", "1.3.1", "1.3.2", "1.4.1", "1.4.2",
91  "1.5.1", "1.5.2", "1.6.1", "1.6.2", "1.7.1", "1.7.2", "1.8.1", "1.8.2",
92  "2.4.1", "2.4.2", "2.5.1", "2.5.2"
93  };
94  for (auto& it : mod) m_PXDModules.push_back(VxdID(it));
95  // set some default size to nu, nv?
96  } else {
97  VXD::SensorInfoBase cellGetInfo = geo.getSensorInfo(m_PXDModules[0]);
98  nu = cellGetInfo.getUCells();
99  nv = cellGetInfo.getVCells();
100  }
101 
102 #ifdef _BELLE2_EPICS
103  m_useEpics |= m_useEpicsRO; // implicit
104  if (m_useEpics) {
105  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
106  }
107 #endif
108 
109  for (VxdID& aPXDModule : m_PXDModules) {
110  TString buff = (std::string)aPXDModule;
111  buff.ReplaceAll(".", "_");
112 #ifdef _BELLE2_EPICS
113  if (m_useEpics) {
114  B2INFO("Connect PVs for " + m_pvPrefix + (std::string)aPXDModule);
115  SEVCHK(ca_create_channel((m_pvPrefix + buff).Data(), NULL, NULL, 10, &mychid_eff[aPXDModule]), "ca_create_channel failure");
116  SEVCHK(ca_create_channel((m_pvPrefix + buff + ".LOW").Data(), NULL, NULL, 10, &mychid_low[aPXDModule]),
117  "ca_create_channel failure");
118  SEVCHK(ca_create_channel((m_pvPrefix + buff + ".LOLO").Data(), NULL, NULL, 10, &mychid_lolo[aPXDModule]),
119  "ca_create_channel failure");
120  }
121 #endif
122  TString histTitle = "PXD Hit Efficiency on Module " + (std::string)aPXDModule + ";Pixel in U;Pixel in V";
123  m_cEffModules[aPXDModule] = new TCanvas((m_histogramDirectoryName + "/c_Eff_").data() + buff);
124  m_hEffModules[aPXDModule] = new TEfficiency("ePXDHitEff_" + buff, histTitle,
125  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
126  }
127 
128  m_cInnerMap = new TCanvas((m_histogramDirectoryName + "/c_InnerMap").data());
129  m_cOuterMap = new TCanvas((m_histogramDirectoryName + "/c_OuterMap").data());
130  m_hInnerMap = new TH2F("hEffInnerMap", "hEffInnerMap", m_u_bins * 8, 0, m_u_bins * 8, m_v_bins * 2, 0, m_v_bins * 2);
131  m_hOuterMap = new TH2F("hEffOuterMap", "hEffOuterMap", m_u_bins * 12, 0, m_u_bins * 12, m_v_bins * 2, 0, m_v_bins * 2);
132 
133  m_hErrorLine = new TH1F("hPXDErrorlimit", "Error Limit", m_PXDModules.size(), 0, m_PXDModules.size());
134  m_hWarnLine = new TH1F("hPXDWarnlimit", "Warn Limit", m_PXDModules.size(), 0, m_PXDModules.size());
135  for (int i = 0; i < (int)m_PXDModules.size(); i++) {
136  m_hErrorLine->SetBinContent(i + 1, m_errorlevel);
137  m_hWarnLine->SetBinContent(i + 1, m_warnlevel);
138  }
139  m_hWarnLine->SetLineColor(kOrange - 3);
140  m_hWarnLine->SetLineWidth(3);
141  m_hWarnLine->SetLineStyle(4);
142  m_hErrorLine->SetLineColor(kRed + 3);
143  m_hErrorLine->SetLineWidth(3);
144  m_hErrorLine->SetLineStyle(7);
145 
146  //One bin for each module in the geometry, one histogram for each layer
147  m_cEffAll = new TCanvas((m_histogramDirectoryName + "/c_EffAll").data());
148  m_hEffAll = new TEfficiency("ePXDHitEffAll", "PXD Integrated Efficiency of each module;PXD Module;",
149  m_PXDModules.size(), 0, m_PXDModules.size());
150  m_hEffAll->SetConfidenceLevel(m_confidence);
151  m_hEffAll->Paint("AP");
152  m_hEffAllLastTotal = m_hEffAll->GetCopyTotalHisto();
153  m_hEffAllLastPassed = m_hEffAll->GetCopyPassedHisto();
154 
155  {
156  auto gr = m_hEffAll->GetPaintedGraph();
157 
158  if (gr) {
159  auto ax = gr->GetXaxis();
160  if (ax) {
161  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
162  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
163  TString ModuleName = (std::string)m_PXDModules[i];
164  ax->SetBinLabel(i + 1, ModuleName);
165  }
166  }
167  }
168  }
169 
170  m_cEffAllUpdate = new TCanvas((m_histogramDirectoryName + "/c_EffAllUp").data());
171  m_hEffAllUpdate = new TEfficiency("ePXDHitEffAllUpdate", "PXD Integral and last-updated Efficiency per module;PXD Module;",
172  m_PXDModules.size(), 0, m_PXDModules.size());
173  m_hEffAllUpdate->SetConfidenceLevel(m_confidence);
174 
175  {
176  auto gr = m_hEffAllUpdate->GetPaintedGraph();
177 
178  if (gr) {
179  auto ax = gr->GetXaxis();
180  if (ax) {
181  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
182  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
183  TString ModuleName = (std::string)m_PXDModules[i];
184  ax->SetBinLabel(i + 1, ModuleName);
185  }
186  }
187  }
188  }
189 
190  m_monObj->addCanvas(m_cEffAll);
191  m_monObj->addCanvas(m_cEffAllUpdate);
192 
193 #ifdef _BELLE2_EPICS
194  if (m_useEpics) {
195  // values per module, see above
196  mychid_status.resize(2);
197  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
198  SEVCHK(ca_create_channel((m_pvPrefix + "Status").data(), NULL, NULL, 10, &mychid_status[0]), "ca_create_channel failure");
199  SEVCHK(ca_create_channel((m_pvPrefix + "Overall").data(), NULL, NULL, 10, &mychid_status[1]), "ca_create_channel failure");
200  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
201  }
202 #endif
203  B2DEBUG(1, "DQMHistAnalysisPXDEff: initialized.");
204 }
205 
206 
207 void DQMHistAnalysisPXDEffModule::beginRun()
208 {
209  B2DEBUG(1, "DQMHistAnalysisPXDEff: beginRun called.");
210 
211  m_cEffAll->Clear();
212  m_cEffAllUpdate->Clear();
213 
214  // no way to reset TEfficiency, do it bin by bin
215  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
216  int j = i + 1;
217  m_hEffAll->SetPassedEvents(j, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
218  m_hEffAll->SetTotalEvents(j, 0);
219  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
220  m_hEffAllUpdate->SetTotalEvents(j, 0);
221 
222  m_warnlevelmod[m_PXDModules[i]] = m_warnlevel;
223  m_errorlevelmod[m_PXDModules[i]] = m_errorlevel;
224 
225 #ifdef _BELLE2_EPICS
226  if (m_useEpics) {
227  // get warn and error limit
228  // as the same array as above, we assume chid exists
229  dbr_double_t tPvData;
230  auto r = ca_get(DBR_DOUBLE, mychid_lolo[m_PXDModules[i]], &tPvData);
231  if (r == ECA_NORMAL) r = ca_pend_io(5.0);// value is only updated here and valid if ECA_NORMAL
232  if (r == ECA_NORMAL) {
233  if (!std::isnan(tPvData)
234  && tPvData > 0.0) {
235  m_hErrorLine->SetBinContent(i + 1, tPvData);
236  if (m_perModuleAlarm) m_errorlevelmod[m_PXDModules[i]] = tPvData;
237  }
238  } else {
239  SEVCHK(r, "ca_get or ca_pend_io failure");
240  }
241  r = ca_get(DBR_DOUBLE, mychid_low[m_PXDModules[i]], &tPvData);
242  if (r == ECA_NORMAL) r = ca_pend_io(5.0);// value is only updated here and valid if ECA_NORMAL
243  if (r == ECA_NORMAL) {
244  if (!std::isnan(tPvData)
245  && tPvData > 0.0) {
246  m_hWarnLine->SetBinContent(i + 1, tPvData);
247  if (m_perModuleAlarm) m_warnlevelmod[m_PXDModules[i]] = tPvData;
248  }
249  } else {
250  SEVCHK(r, "ca_get or ca_pend_io failure");
251  }
252  }
253 #endif
254 
255  }
256  // Thus histo will contain old content until first update
257  m_hEffAllLastTotal->Reset();
258  m_hEffAllLastPassed->Reset();
259 
260  for (auto single_cmap : m_cEffModules) {
261  if (single_cmap.second) single_cmap.second->Clear();
262  }
263 
264 }
265 
266 
267 void DQMHistAnalysisPXDEffModule::event()
268 {
269  //Count how many of each type of histogram there are for the averaging
270  //std::map<std::string, int> typeCounter;
271 
272  m_hInnerMap->Reset();
273  m_hOuterMap->Reset();
274 
275  for (unsigned int i = 1; i <= m_PXDModules.size(); i++) {
276  VxdID& aPXDModule = m_PXDModules[i - 1];
277 
278  TString buff = (std::string)aPXDModule;
279  buff.ReplaceAll(".", "_");
280 
281  TH1* Hits, *Matches;
282  TString locationHits = "track_hits_" + buff;
283  if (m_histogramDirectoryName != "") {
284  locationHits = m_histogramDirectoryName + "/" + locationHits;
285  }
286  Hits = (TH1*)findHist(locationHits.Data());
287  TString locationMatches = "matched_cluster_" + buff;
288  if (m_histogramDirectoryName != "") {
289  locationMatches = m_histogramDirectoryName + "/" + locationMatches;
290  }
291  Matches = (TH1*)findHist(locationMatches.Data());
292 
293  // Finding only one of them should only happen in very strange situations...
294  if (Hits == nullptr || Matches == nullptr) {
295  B2ERROR("Missing histogram for sensor " << aPXDModule);
296  } else {
297  if (m_cEffModules[aPXDModule] && m_hEffModules[aPXDModule]) {// this check creates them with a nullptr ..bad
298  m_hEffModules[aPXDModule]->SetTotalHistogram(*Hits, "f");
299  m_hEffModules[aPXDModule]->SetPassedHistogram(*Matches, "f");
300 
301  m_cEffModules[aPXDModule]->cd();
302  m_hEffModules[aPXDModule]->Paint("colz"); // not Draw, enforce to create GetPaintedHistogram?
303  m_cEffModules[aPXDModule]->Modified();
304  m_cEffModules[aPXDModule]->Update();
305 
306  auto h = m_hEffModules[aPXDModule]->GetPaintedHistogram();
307  int s = (2 - aPXDModule.getSensorNumber()) * m_v_bins;
308  int l = (aPXDModule.getLadderNumber() - 1) * m_u_bins;
309  if (m_hInnerMap && aPXDModule.getLayerNumber() == 1) {
310  for (int u = 0; u < m_u_bins; u++) {
311  for (int v = 0; v < m_v_bins; v++) {
312  auto b = h->GetBin(u + 1, v + 1);
313  m_hInnerMap->Fill(u + l, v + s, h->GetBinContent(b));
314  }
315  }
316  }
317  if (m_hOuterMap && aPXDModule.getLayerNumber() == 2) {
318  for (int u = 0; u < m_u_bins; u++) {
319  for (int v = 0; v < m_v_bins; v++) {
320  auto b = h->GetBin(u + 1, v + 1);
321  m_hOuterMap->Fill(u + l, v + s, h->GetBinContent(b));
322  }
323  }
324  }
325  }
326  }
327  }// Single-Module histos + 2d overview finished
328 
329  m_cInnerMap->cd();
330  m_hInnerMap->Draw("colz");
331  m_cInnerMap->Modified();
332  m_cInnerMap->Update();
333  m_cOuterMap->cd();
334  m_hOuterMap->Draw("colz");
335  m_cOuterMap->Modified();
336  m_cOuterMap->Update();
337 
338  // Change: We now use one histogram for hits and matches to make
339  // sure that we have an atomic update which is otherwise not
340  // guaranteed by DQM framework
341  TString locationHits = "PXD_Eff_combined";
342  if (m_histogramDirectoryName != "") {
343  locationHits = m_histogramDirectoryName + "/" + locationHits;
344  }
345  TH1* Combined = (TH1*)findHist(locationHits.Data());
346 
347  double stat_data = 0;
348  bool error_flag = false;
349  bool warn_flag = false;
350  double all = 0.0;
351 
352  double imatch = 0.0, ihit = 0.0;
353  int ieff = 0;
354 
355  std::map <VxdID, bool> updated{}; // init to false, keep track of updated histograms
356  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
357  int j = i + 1;
358 
359  if (Combined == nullptr) {
360  m_hEffAll->SetPassedEvents(j, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
361  m_hEffAll->SetTotalEvents(j, 0);
362  } else {
363  VxdID& aModule = m_PXDModules[i];
364  double nmatch = Combined->GetBinContent(i * 2 + 2);
365  double nhit = Combined->GetBinContent(i * 2 + 1);
366  if (nmatch > 10 && nhit > 10) { // could be zero, too
367  imatch += nmatch;
368  ihit += nhit;
369  ieff++; // only count in modules working
370  double var_e = nmatch / nhit; // can never be zero
371  // keep workaround code for re-use
372 // if (j == 6) continue; // workaround for 1.3.2 module
373  m_monObj->setVariable(Form("efficiency_%d_%d_%d", aModule.getLayerNumber(), aModule.getLadderNumber(), aModule.getSensorNumber()),
374  var_e);
375  }
376 
378  all += nhit;
379  m_hEffAll->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
380  m_hEffAll->SetTotalEvents(j, nhit);
381  m_hEffAll->SetPassedEvents(j, nmatch);
382 
383  if (nhit < m_minEntries) {
384  // update the first entries directly (short runs)
385  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
386  m_hEffAllUpdate->SetTotalEvents(j, nhit);
387  m_hEffAllUpdate->SetPassedEvents(j, nmatch);
388  m_hEffAllLastTotal->SetBinContent(j, nhit);
389  m_hEffAllLastPassed->SetBinContent(j, nmatch);
390  updated[aModule] = true;
391  } else if (nhit - m_hEffAllLastTotal->GetBinContent(j) > m_minEntries) {
392  m_hEffAllUpdate->SetPassedEvents(j, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
393  m_hEffAllUpdate->SetTotalEvents(j, nhit - m_hEffAllLastTotal->GetBinContent(j));
394  m_hEffAllUpdate->SetPassedEvents(j, nmatch - m_hEffAllLastPassed->GetBinContent(j));
395  m_hEffAllLastTotal->SetBinContent(j, nhit);
396  m_hEffAllLastPassed->SetBinContent(j, nmatch);
397  updated[aModule] = true;
398  }
399 
400  // keep workaround code for re-use
401 // if (j == 6) continue; // workaround for 1.3.2 module
402 
403  // get the errors and check for limits for each bin seperately ...
404 
405  if (nhit > 50) {
406  error_flag |= (m_hEffAll->GetEfficiency(j) + m_hEffAll->GetEfficiencyErrorUp(j) <
407  m_errorlevelmod[aModule]); // error if upper error value is below limit
408  warn_flag |= (m_hEffAll->GetEfficiency(j) + m_hEffAll->GetEfficiencyErrorUp(j) <
409  m_warnlevelmod[aModule]); // (and not only the actual eff value)
410  if (m_alarmAdhoc) {
411  error_flag |= (m_hEffAllUpdate->GetEfficiency(j) + m_hEffAllUpdate->GetEfficiencyErrorUp(j) <
412  m_errorlevelmod[aModule]); // error if upper error value is below limit
413  warn_flag |= (m_hEffAllUpdate->GetEfficiency(j) + m_hEffAllUpdate->GetEfficiencyErrorUp(j) <
414  m_warnlevelmod[aModule]); // (and not only the actual eff value)
415  }
416  }
417  }
418  }
419 
420  {
421  m_cEffAll->cd();
422  m_cEffAll->cd(0);
423  m_hEffAll->Paint("AP");
424  m_cEffAll->Clear();
425  m_cEffAll->cd(0);
426 
427  auto gr = m_hEffAll->GetPaintedGraph();
428  if (gr) {
429  double scale_min = 1.0;
430  for (int i = 0; i < gr->GetN(); i++) {
431  gr->SetPointEXhigh(i, 0.);
432  gr->SetPointEXlow(i, 0.);
433  // this has to be done first, as it will recalc Min/Max and destroy axis
434  Double_t x, y;
435  gr->GetPoint(i, x, y);
436  gr->SetPoint(i, x - 0.01, y); // workaround for jsroot bug (fixed upstream)
437  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
438  // keep workaround code for re-use
439 // if (i != 5) // exclude 1.3.2
440  if (scale_min > val) scale_min = val;
441  }
442  if (scale_min == 1.0) scale_min = 0.0;
443  if (scale_min > 0.9) scale_min = 0.9;
444  gr->SetMinimum(0);
445  gr->SetMaximum(m_PXDModules.size());
446  auto ay = gr->GetYaxis();
447  if (ay) ay->SetRangeUser(scale_min, 1.0);
448  auto ax = gr->GetXaxis();
449  if (ax) {
450  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
451  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
452  TString ModuleName = (std::string)m_PXDModules[i];
453  ax->SetBinLabel(i + 1, ModuleName);
454  }
455  }
456 
457  gr->SetLineColor(4);
458  gr->SetLineWidth(2);
459  gr->SetMarkerStyle(8);
460 
461  gr->Draw("AP");
462 
463  // keep workaround code for re-use
464 // auto tt = new TLatex(5.5, scale_min, " 1.3.2 Module is excluded, please ignore");
465 // tt->SetTextAngle(90);// Rotated
466 // tt->SetTextAlign(12);// Centered
467 // tt->Draw();
468 
469  if (all < 100.) {
470  m_cEffAll->Pad()->SetFillColor(kGray);// Magenta or Gray
471  } else {
472  if (error_flag) {
473  m_cEffAll->Pad()->SetFillColor(kRed);// Red
474  } else if (warn_flag) {
475  m_cEffAll->Pad()->SetFillColor(kYellow);// Yellow
476  } else {
477  m_cEffAll->Pad()->SetFillColor(kGreen);// Green
478  // m_cEffAll->Pad()->SetFillColor(kWhite);// White
479  }
480  }
481 
482  m_cEffAll->Pad()->SetFrameFillColor(kWhite - 1); // White
483  m_cEffAll->Pad()->SetFrameFillStyle(1001);// White
484  m_cEffAll->Pad()->Modified();
485  m_cEffAll->Pad()->Update();
486  m_hWarnLine->Draw("same,hist");
487  m_hErrorLine->Draw("same,hist");
488  }
489 
490  m_cEffAll->Modified();
491  m_cEffAll->Update();
492  }
493 
494  {
495  m_cEffAllUpdate->cd();
496  m_hEffAllUpdate->Paint("AP");
497  m_cEffAllUpdate->Clear();
498  m_cEffAllUpdate->cd(0);
499 
500  auto gr = m_hEffAllUpdate->GetPaintedGraph();
501  auto gr3 = (TGraphAsymmErrors*) m_hEffAll->GetPaintedGraph()->Clone();
502  if (gr3) {
503  for (int i = 0; i < gr3->GetN(); i++) {
504  Double_t x, y;
505  gr3->GetPoint(i, x, y);
506  gr3->SetPoint(i, x + 0.2, y);
507  }
508  }
509 
510  double scale_min = 1.0;
511  if (gr) {
512  for (int i = 0; i < gr->GetN(); i++) {
513  gr->SetPointEXhigh(i, 0.);
514  gr->SetPointEXlow(i, 0.);
515  // this has to be done first, as it will recalc Min/Max and destroy axis
516  Double_t x, y;
517  gr->GetPoint(i, x, y);
518  gr->SetPoint(i, x - 0.2, y); // shift a bit if in same plot
519  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
520  // keep workaround code for re-use
521 // if (i != 5) { // exclude 1.3.2
522  if (scale_min > val) scale_min = val;
523  }
524  if (scale_min == 1.0) scale_min = 0.0;
525  if (scale_min > 0.9) scale_min = 0.9;
526  gr->SetMinimum(0);
527  gr->SetMaximum(m_PXDModules.size());
528  auto ay = gr->GetYaxis();
529  if (ay) ay->SetRangeUser(scale_min, 1.0);
530  auto ax = gr->GetXaxis();
531  if (ax) {
532  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
533  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
534  TString ModuleName = (std::string)m_PXDModules[i];
535  ax->SetBinLabel(i + 1, ModuleName);
536  }
537  }
538 #ifdef _BELLE2_EPICS
539  if (m_useEpics && !m_useEpicsRO) {
540  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
541  if (updated[m_PXDModules[i]]) {
542  Double_t x, y;// we assume that double and Double_t are same!
543  gr->GetPoint(i, x, y);
544  auto& my = mychid_eff[m_PXDModules[i]];// as the same array as above, we assume it exists
545  // we should only write if it was updated!
546  SEVCHK(ca_put(DBR_DOUBLE, my, (void*)&y), "ca_set failure");
547  }
548  }
549  }
550 #endif
551  gr->SetLineColor(kBlack);
552  gr->SetLineWidth(3);
553  gr->SetMarkerStyle(33);
554  } else scale_min = 0.0;
555  if (gr) gr->Draw("AP");
556  if (gr3) gr3->Draw("P");
557  // keep workaround code for re-use
558 // auto tt = new TLatex(5.5, scale_min, " 1.3.2 Module is excluded, please ignore");
559 // tt->SetTextSize(0.035);
560 // tt->SetTextAngle(90);// Rotated
561 // tt->SetTextAlign(12);// Centered
562 // tt->Draw();
563 
564  if (all < 100.) {
565  m_cEffAllUpdate->Pad()->SetFillColor(kGray);// Magenta or Gray
566  stat_data = 0.;
567  } else {
568  if (error_flag) {
569  m_cEffAllUpdate->Pad()->SetFillColor(kRed);// Red
570  stat_data = 4.;
571  } else if (warn_flag) {
572  m_cEffAllUpdate->Pad()->SetFillColor(kYellow);// Yellow
573  stat_data = 3.;
574  } else {
575  m_cEffAllUpdate->Pad()->SetFillColor(kGreen);// Green
576  stat_data = 2.;
578  // m_cEffAllUpdate->Pad()->SetFillColor(kWhite);// White
579  }
580  }
581  m_cEffAllUpdate->Pad()->SetFrameFillColor(kWhite - 1); // White
582  m_cEffAllUpdate->Pad()->SetFrameFillStyle(1001);// White
583  m_cEffAllUpdate->Pad()->Modified();
584  m_cEffAllUpdate->Pad()->Update();
585  m_hWarnLine->Draw("same,hist");
586  m_hErrorLine->Draw("same,hist");
587  }
588  m_cEffAllUpdate->Modified();
589  m_cEffAllUpdate->Update();
590 
591 
592  double var_efficiency = ihit > 0 ? imatch / ihit : 0.0;
593  m_monObj->setVariable("efficiency", var_efficiency);
594  m_monObj->setVariable("nmodules", ieff);
595 
596 #ifdef _BELLE2_EPICS
597  if (m_useEpics && !m_useEpicsRO) {
598  SEVCHK(ca_put(DBR_DOUBLE, mychid_status[0], (void*)&stat_data), "ca_set failure");
599  // only update if statistics is reasonable, we dont want "0" drops between runs!
600  if (stat_data != 0) SEVCHK(ca_put(DBR_DOUBLE, mychid_status[1], (void*)&var_efficiency), "ca_set failure");
601  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
602  }
603 #endif
604 }
605 
606 void DQMHistAnalysisPXDEffModule::terminate()
607 {
608  B2DEBUG(1, "DQMHistAnalysisPXDEff: terminate called");
609 #ifdef _BELLE2_EPICS
610  if (m_useEpics) {
611  for (auto& m : mychid_status) SEVCHK(ca_clear_channel(m), "ca_clear_channel failure");
612  for (auto& m : mychid_eff) SEVCHK(ca_clear_channel(m.second), "ca_clear_channel failure");
613  for (auto& m : mychid_low) SEVCHK(ca_clear_channel(m.second), "ca_clear_channel failure");
614  for (auto& m : mychid_lolo) SEVCHK(ca_clear_channel(m.second), "ca_clear_channel failure");
615  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
616  }
617 #endif
618 }
619 
The base class for the histogram analysis module.
DQM Histogram Analysis for PXD Efficiency.
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:58
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
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
#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.