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