Belle II Software  release-08-01-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.99);
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, "the list of excluded modules, indices 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_eEffModules[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_nrxbins = m_PXDModules.size() + 3; // Modules + L1 + L2 + All
115  m_hErrorLine = new TH1F("hPXDErrorlimit", "Error Limit", m_nrxbins, 0, m_nrxbins);
116  m_hWarnLine = new TH1F("hPXDWarnlimit", "Warn Limit", m_nrxbins, 0, m_nrxbins);
117  for (int i = 0; i < (int)m_nrxbins; i++) {
118  m_hErrorLine->SetBinContent(i + 1, m_errorlevel);
119  m_hWarnLine->SetBinContent(i + 1, m_warnlevel);
120  }
121  m_hWarnLine->SetLineColor(kOrange - 3);
122  m_hWarnLine->SetLineWidth(3);
123  m_hWarnLine->SetLineStyle(4);
124  m_hErrorLine->SetLineColor(kRed + 3);
125  m_hErrorLine->SetLineWidth(3);
126  m_hErrorLine->SetLineStyle(7);
127 
128  //One bin for each module in the geometry, one histogram for each layer
129  m_cEffAll = new TCanvas((m_histogramDirectoryName + "/c_EffAll").data());
130  m_eEffAll = new TEfficiency("ePXDHitEffAll", "PXD Integrated Efficiency of each module;PXD Module;", m_nrxbins, 0, m_nrxbins);
131  m_eEffAll->SetConfidenceLevel(m_confidence);
132  m_eEffAll->Paint("AP");
133  m_hEffAllLastTotal = m_eEffAll->GetCopyTotalHisto();
134  m_hEffAllLastPassed = m_eEffAll->GetCopyPassedHisto();
135 
136  setLabels(m_eEffAll->GetPaintedGraph());
137 
138  m_cEffAllUpdate = new TCanvas((m_histogramDirectoryName + "/c_EffAllUp").data());
139  m_eEffAllUpdate = new TEfficiency("ePXDHitEffAllUpdate", "PXD Integral and last-updated Efficiency per module;PXD Module;",
140  m_nrxbins, 0, m_nrxbins);
141  m_eEffAllUpdate->SetConfidenceLevel(m_confidence);
142 
143  m_eEffAllUpdate->Paint("AP");
144  setLabels(m_eEffAllUpdate->GetPaintedGraph());
145 
148 
149  registerEpicsPV("PXD:Eff:Status", "Status");
150  registerEpicsPV("PXD:Eff:Overall", "All");
151  registerEpicsPV("PXD:Eff:L1", "L1");
152  registerEpicsPV("PXD:Eff:L2", "L2");
153  B2DEBUG(1, "DQMHistAnalysisPXDEff: initialized.");
154 }
155 
156 
158 {
159  B2DEBUG(1, "DQMHistAnalysisPXDEff: beginRun called.");
160 
161  // Clear all used canvases
162  m_cEffAll->Clear();
163  m_cEffAllUpdate->Clear();
164  m_cInnerMap->Clear();
165  m_cOuterMap->Clear();
166  for (auto single_cmap : m_cEffModules) {
167  if (single_cmap.second) single_cmap.second->Clear();
168  }
169 
170  // The 2d Efficiency maps (m_eEffModules[]) per module are not cleared, but re-created each update
171  // also they are only drawn to Canvas on update, thus no clear is needed here
172 
173  // Reset TEfficiency and get (new) alarm limits from PVs
174  // no way to reset TEfficiency, do it bin by bin
175  for (int i = 0; i < m_nrxbins; i++) {
176  int bin = i + 1;
177  m_eEffAll->SetPassedEvents(bin, 0); // order, otherwise it might happen that SetTotalEvents is NOT filling the value!
178  m_eEffAll->SetTotalEvents(bin, 0);
179  m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
180  m_eEffAllUpdate->SetTotalEvents(bin, 0);
181 
182  if (i < int(m_PXDModules.size())) { // only for modules
185 
186  // get warn and error limit
187  // as the same array as above, we assume chid exists
188  double dummy, loerr = 0, lowarn = 0;
189  if (requestLimitsFromEpicsPVs((std::string)m_PXDModules[i], loerr, lowarn, dummy, dummy)) {
190  m_hErrorLine->SetBinContent(bin, loerr);
192  m_hWarnLine->SetBinContent(bin, lowarn);
193  if (m_perModuleAlarm) m_warnlevelmod[m_PXDModules[i]] = lowarn;
194  }
195  }
196  }
197  {
198  double dummy, loerr = 0, lowarn = 0;
199  m_warnlevelmod["L1"] = m_warnlevel;
201  if (requestLimitsFromEpicsPVs("L1", loerr, lowarn, dummy, dummy)) {
202  m_hErrorLine->SetBinContent(m_PXDModules.size() + 1, loerr);
203  if (m_perModuleAlarm) m_errorlevelmod["L1"] = loerr;
204  m_hWarnLine->SetBinContent(m_PXDModules.size() + 1, lowarn);
205  if (m_perModuleAlarm) m_warnlevelmod["L1"] = lowarn;
206  }
207  m_warnlevelmod["L2"] = m_warnlevel;
209  if (requestLimitsFromEpicsPVs("L2", loerr, lowarn, dummy, dummy)) {
210  m_hErrorLine->SetBinContent(m_PXDModules.size() + 2, loerr);
211  if (m_perModuleAlarm) m_errorlevelmod["L2"] = loerr;
212  m_hWarnLine->SetBinContent(m_PXDModules.size() + 2, lowarn);
213  if (m_perModuleAlarm) m_warnlevelmod["L2"] = lowarn;
214  }
215  m_warnlevelmod["All"] = m_warnlevel;
216  m_errorlevelmod["All"] = m_errorlevel;
217  if (requestLimitsFromEpicsPVs("All", loerr, lowarn, dummy, dummy)) {
218  m_hErrorLine->SetBinContent(m_PXDModules.size() + 3, loerr);
219  if (m_perModuleAlarm) m_errorlevelmod["All"] = loerr;
220  m_hWarnLine->SetBinContent(m_PXDModules.size() + 3, lowarn);
221  if (m_perModuleAlarm) m_warnlevelmod["All"] = lowarn;
222  }
223  }
224 
225  // Clear all remaining Histograms (e.g. for our private delta histogramming)
226  m_hEffAllLastTotal->Reset();
227  m_hEffAllLastPassed->Reset();
228  m_hInnerMap->Reset();
229  m_hOuterMap->Reset();
230 }
231 
232 bool DQMHistAnalysisPXDEffModule::updateEffBins(int bin, int nhit, int nmatch, int minentries)
233 {
234  m_eEffAll->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
235  m_eEffAll->SetTotalEvents(bin, nhit);
236  m_eEffAll->SetPassedEvents(bin, nmatch);
237 
238  if (nhit < minentries) {
239  // update the first entries directly (short runs)
240  m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
241  m_eEffAllUpdate->SetTotalEvents(bin, nhit);
242  m_eEffAllUpdate->SetPassedEvents(bin, nmatch);
243  m_hEffAllLastTotal->SetBinContent(bin, nhit);
244  m_hEffAllLastPassed->SetBinContent(bin, nmatch);
245  return true;
246  } else if (nhit - m_hEffAllLastTotal->GetBinContent(bin) > minentries) {
247  m_eEffAllUpdate->SetPassedEvents(bin, 0); // otherwise it might happen that SetTotalEvents is NOT filling the value!
248  m_eEffAllUpdate->SetTotalEvents(bin, nhit - m_hEffAllLastTotal->GetBinContent(bin));
249  m_eEffAllUpdate->SetPassedEvents(bin, nmatch - m_hEffAllLastPassed->GetBinContent(bin));
250  m_hEffAllLastTotal->SetBinContent(bin, nhit);
251  m_hEffAllLastPassed->SetBinContent(bin, nmatch);
252  return true;
253  }// else
254  return false;
255 }
256 
257 bool DQMHistAnalysisPXDEffModule::check_warn_level(int bin, std::string name)
258 {
259  bool warn_flag = (m_eEffAll->GetEfficiency(bin) + m_eEffAll->GetEfficiencyErrorUp(bin) <
260  m_warnlevelmod[name]); // (and not only the actual eff value)
261  if (m_alarmAdhoc) {
262  warn_flag |= (m_eEffAllUpdate->GetEfficiency(bin) + m_eEffAllUpdate->GetEfficiencyErrorUp(bin) <
263  m_warnlevelmod[name]); // (and not only the actual eff value)
264  }
265  return warn_flag;
266 }
267 
268 bool DQMHistAnalysisPXDEffModule::check_error_level(int bin, std::string name)
269 {
270  bool error_flag = (m_eEffAll->GetEfficiency(bin) + m_eEffAll->GetEfficiencyErrorUp(bin) <
271  m_errorlevelmod[name]); // error if upper error value is below limit
272  if (m_alarmAdhoc) {
273  error_flag |= (m_eEffAllUpdate->GetEfficiency(bin) + m_eEffAllUpdate->GetEfficiencyErrorUp(bin) <
274  m_errorlevelmod[name]); // error if upper error value is below limit
275  }
276  return error_flag;
277 }
278 
279 void DQMHistAnalysisPXDEffModule::setLabels(TGraphAsymmErrors* gr)
280 {
281  if (gr) {
282  auto ax = gr->GetXaxis();
283  if (ax) {
284  ax->Set(m_nrxbins, 0, m_nrxbins);
285  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
286  TString ModuleName = (std::string)m_PXDModules[i];
287  ax->SetBinLabel(i + 1, ModuleName);
288  }
289  ax->SetBinLabel(m_PXDModules.size() + 1, "L1");
290  ax->SetBinLabel(m_PXDModules.size() + 2, "L2");
291  ax->SetBinLabel(m_PXDModules.size() + 3, "All");
292  }
293  }
294 }
295 
297 {
298  {
299  // First create some 2d overview of efficiency for all modules
300  // This is not taken into account for efficiency calculation as
301  // there may be update glitches dues to seperate histograms
302  // The histograms
303  bool updateinner = false, updateouter = false;
304  for (auto aPXDModule : m_PXDModules) {
305  auto buff = (std::string)aPXDModule;
306  replace(buff.begin(), buff.end(), '.', '_');
307 
308  std::string locationHits = "track_hits_" + buff;
309  if (m_histogramDirectoryName != "") {
310  locationHits = m_histogramDirectoryName + "/" + locationHits;
311  }
312  std::string locationMatches = "matched_cluster_" + buff;
313  if (m_histogramDirectoryName != "") {
314  locationMatches = m_histogramDirectoryName + "/" + locationMatches;
315  }
316 
317  auto Hits = (TH1*)findHist(locationHits, true);// check if updated
318  auto Matches = (TH1*)findHist(locationMatches, true);// check if updated
319 
320  if (Hits == nullptr && Matches == nullptr) continue; // none updated
321 
322  if (Hits == nullptr) Hits = (TH1*)findHist(locationHits); // actually, this should not happen ...
323  if (Matches == nullptr) Matches = (TH1*)findHist(locationMatches); // ... as updates should coincide
324 
325  // Finding only one of them should only happen in very strange situations... still better check
326  if (Hits && Matches) {
327  if (m_cEffModules[aPXDModule] && m_eEffModules[aPXDModule]) {// this check creates them with a nullptr ..bad
328  m_eEffModules[aPXDModule]->SetTotalHistogram(*Hits, "f");
329  m_eEffModules[aPXDModule]->SetPassedHistogram(*Matches, "f");
330 
331  m_cEffModules[aPXDModule]->cd();
332  m_eEffModules[aPXDModule]->Paint("colz"); // not Draw, enforce to create GetPaintedHistogram?
333  m_eEffModules[aPXDModule]->Draw("colz"); // but Draw needed to export Canvas!
334  m_cEffModules[aPXDModule]->Modified();
335  m_cEffModules[aPXDModule]->Update();
336  UpdateCanvas(m_cEffModules[aPXDModule]);
337 
338  auto h = m_eEffModules[aPXDModule]->GetPaintedHistogram();
339  int s = (2 - aPXDModule.getSensorNumber()) * m_v_bins;
340  int l = (aPXDModule.getLadderNumber() - 1) * m_u_bins;
341  if (m_hInnerMap && aPXDModule.getLayerNumber() == 1) {
342  updateinner = true;
343  for (int u = 0; u < m_u_bins; u++) {
344  for (int v = 0; v < m_v_bins; v++) {
345  auto b = h->GetBin(u + 1, v + 1);
346  m_hInnerMap->Fill(u + l, v + s, h->GetBinContent(b));
347  }
348  }
349  }
350  if (m_hOuterMap && aPXDModule.getLayerNumber() == 2) {
351  updateouter = true;
352  for (int u = 0; u < m_u_bins; u++) {
353  for (int v = 0; v < m_v_bins; v++) {
354  auto b = h->GetBin(u + 1, v + 1);
355  m_hOuterMap->Fill(u + l, v + s, h->GetBinContent(b));
356  }
357  }
358  }
359  }
360  } else {
361  B2WARNING("only one plot upd " << aPXDModule);
362  }
363  }
364  // Single-Module histos + 2d overview finished. now draw overviews
365  if (updateinner) {
366  m_cInnerMap->cd();
367  if (m_hInnerMap) m_hInnerMap->Draw("colz");
368  m_cInnerMap->Modified();
369  m_cInnerMap->Update();
371  }
372  if (updateouter) {
373  m_cOuterMap->cd();
374  if (m_hOuterMap) m_hOuterMap->Draw("colz");
375  m_cOuterMap->Modified();
376  m_cOuterMap->Update();
378  }
379  // 3d overview done
380  }
381 
382 
383 // Now, calculate and update efficiency.
384 // Only histogram is used for numerator AND denominatorto to have an atomic update.
385 // (avoid possible update glitches from daq/dqm framework side)
386 // The bins per module are read out and filled into an TEfficiency as total and passed events into bin
387 // Summaries for L1, L2, Overall are added, too
388  auto Combined = (TH1*)findHist(m_histogramDirectoryName, "PXD_Eff_combined", true);// only if updated
389 
390  if (Combined) {
391  // only if histogram was changed
392 
393  EStatus stat_data = c_StatusTooFew;
394  bool error_flag = false;
395  bool warn_flag = false;
396  double all = 0.0;
397 
398  double imatch = 0.0, ihit = 0.0;
399  double imatchL1 = 0.0, ihitL1 = 0.0;
400  double imatchL2 = 0.0, ihitL2 = 0.0;
401  int ieff = 0; // count number of modules with useful stytistics
402 
403  std::map <VxdID, bool> updated{}; // init to false, keep track of updated histograms
404  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
405  // workaround for excluded module
406  if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
407  // excluded modules are not counted at all!
408  int bin = i + 1; // bin nr is index +1
409 
410  VxdID& aModule = m_PXDModules[i];
411  double nmatch = Combined->GetBinContent(i * 2 + 2);
412  double nhit = Combined->GetBinContent(i * 2 + 1);
413 
414  imatch += nmatch;
415  ihit += nhit;
416  // check layer
417  if (i >= 16) {
418  imatchL2 += nmatch;
419  ihitL2 += nhit;
420  } else {
421  imatchL1 += nmatch;
422  ihitL1 += nhit;
423  }
424 
425  if (nhit >= m_minEntries) { // dont update if there is nothing to calculate
426  ieff++; // only count in modules with significant stat
427  double var_e = nmatch / nhit; // can never be zero
428  m_monObj->setVariable(Form("efficiency_%d_%d_%d", aModule.getLayerNumber(), aModule.getLadderNumber(), aModule.getSensorNumber()),
429  var_e);
430  }
431 
433  all += nhit;
434 
435  updated[aModule] = updateEffBins(bin, nhit, nmatch, m_minEntries);
436 
437  // workaround for excluded module
438  if (std::find(m_excluded.begin(), m_excluded.end(), i) != m_excluded.end()) continue;
439 
440  // get the errors and check for limits for each bin seperately ...
441 
442  if (nhit >= m_minEntries) {
443  error_flag |= check_error_level(bin, aModule);
444  warn_flag |= check_warn_level(bin, aModule);
445  }
446  }
447 
448  updateEffBins(m_PXDModules.size() + 1, ihitL1, imatchL1, m_minEntries * 8);
449  if (ihitL1 >= m_minEntries) {
450  error_flag |= check_error_level(m_PXDModules.size() + 1, "L1");
451  warn_flag |= check_warn_level(m_PXDModules.size() + 1, "L1");
452  }
453  updateEffBins(m_PXDModules.size() + 2, ihitL2, imatchL2, m_minEntries * 12);
454  if (ihitL2 >= m_minEntries) {
455  error_flag |= check_error_level(m_PXDModules.size() + 2, "L2");
456  warn_flag |= check_warn_level(m_PXDModules.size() + 2, "L2");
457  }
458  updateEffBins(m_PXDModules.size() + 3, ihit, imatch, m_minEntries * 20);
459  if (ihit >= m_minEntries) {
460  error_flag |= check_error_level(m_PXDModules.size() + 3, "All");
461  warn_flag |= check_warn_level(m_PXDModules.size() + 3, "All");
462  }
463 
464  {
465  m_cEffAll->cd();
466  m_cEffAll->cd(0);
467  m_eEffAll->Paint("AP");
468  m_cEffAll->Clear();
469  m_cEffAll->cd(0);
470 
471  auto gr = m_eEffAll->GetPaintedGraph();
472  if (gr) {
473  double scale_min = 1.0;
474  for (int i = 0; i < gr->GetN(); i++) {
475  gr->SetPointEXhigh(i, 0.);
476  gr->SetPointEXlow(i, 0.);
477  // this has to be done first, as it will recalc Min/Max and destroy axis
478  Double_t x, y;
479  gr->GetPoint(i, x, y);
480  gr->SetPoint(i, x - 0.01, y); // workaround for jsroot bug (fixed upstream)
481  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
482  if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
483  // scale update only for included module
484  if (scale_min > val) scale_min = val;
485  }
486  }
487  if (scale_min == 1.0) scale_min = 0.0;
488  if (scale_min > 0.9) scale_min = 0.9;
489  auto ay = gr->GetYaxis();
490  if (ay) ay->SetRangeUser(scale_min, 1.0);
491  setLabels(gr);
492 
493  gr->SetLineColor(4);
494  gr->SetLineWidth(2);
495  gr->SetMarkerStyle(8);
496 
497  gr->Draw("AP");
498 
499  for (auto& it : m_excluded) {
500  auto tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
501  tt->SetTextSize(0.035);
502  tt->SetTextAngle(90);// Rotated
503  tt->SetTextAlign(12);// Centered
504  tt->Draw();
505  }
506 
507 
508  EStatus all_stat = makeStatus(all >= m_minEntries, warn_flag, error_flag);
509  colorizeCanvas(m_cEffAll, all_stat);
510 
511  m_hWarnLine->Draw("same,hist");
512  m_hErrorLine->Draw("same,hist");
513  }
514 
515  UpdateCanvas(m_cEffAll->GetName());
516  m_cEffAll->Modified();
517  m_cEffAll->Update();
518  }
519 
520  {
521  m_cEffAllUpdate->cd();
522  m_eEffAllUpdate->Paint("AP");
523  m_cEffAllUpdate->Clear();
524  m_cEffAllUpdate->cd(0);
525 
526  auto gr = m_eEffAllUpdate->GetPaintedGraph();
527  auto gr3 = (TGraphAsymmErrors*) m_eEffAll->GetPaintedGraph()->Clone();
528  if (gr3) {
529  for (int i = 0; i < gr3->GetN(); i++) {
530  Double_t x, y;
531  gr3->GetPoint(i, x, y);
532  gr3->SetPoint(i, x + 0.2, y);
533  }
534  }
535 
536  double scale_min = 1.0;
537  if (gr) {
538  for (int i = 0; i < gr->GetN(); i++) {
539  gr->SetPointEXhigh(i, 0.);
540  gr->SetPointEXlow(i, 0.);
541  // this has to be done first, as it will recalc Min/Max and destroy axis
542  Double_t x, y;
543  gr->GetPoint(i, x, y);
544  gr->SetPoint(i, x - 0.2, y); // shift a bit if in same plot
545  auto val = y - gr->GetErrorYlow(i); // Error is relative to value
546  if (std::find(m_excluded.begin(), m_excluded.end(), i) == m_excluded.end()) {
547  // skip scale update only for included modules
548  if (scale_min > val) scale_min = val;
549  }
550  }
551  if (scale_min == 1.0) scale_min = 0.0;
552  if (scale_min > 0.9) scale_min = 0.9;
553  auto ay = gr->GetYaxis();
554  if (ay) ay->SetRangeUser(scale_min, 1.0);
555  setLabels(gr);
556 
557  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
558  if (updated[m_PXDModules[i]]) {
559  // we should only write if it was updated!
560  Double_t x, y;// we assume that double and Double_t are same!
561  gr->GetPoint(i, x, y);
562  setEpicsPV((std::string)m_PXDModules[i], y);
563  }
564  }
565 
566  gr->SetLineColor(kBlack);
567  gr->SetLineWidth(3);
568  gr->SetMarkerStyle(33);
569  gr->Draw("AP");
570  } else scale_min = 0.0;
571  if (gr3) gr3->Draw("P"); // both in one plot
572 
573  for (auto& it : m_excluded) {
574  auto tt = new TLatex(it + 0.5, scale_min, (" " + std::string(m_PXDModules[it]) + " Module is excluded, please ignore").c_str());
575  tt->SetTextSize(0.035);
576  tt->SetTextAngle(90);// Rotated
577  tt->SetTextAlign(12);// Centered
578  tt->Draw();
579  }
580 
581  stat_data = makeStatus(all >= m_minEntries, warn_flag, error_flag);
582  colorizeCanvas(m_cEffAllUpdate, stat_data);
583 
584  m_hWarnLine->Draw("same,hist");
585  m_hErrorLine->Draw("same,hist");
586  }
587  UpdateCanvas(m_cEffAllUpdate->GetName());
588  m_cEffAllUpdate->Modified();
589  m_cEffAllUpdate->Update();
590 
591  double var_efficiency = ihit > 0 ? imatch / ihit : 0.0;
592  double var_efficiencyL1 = ihitL1 > 0 ? imatchL1 / ihitL1 : 0.0;
593  double var_efficiencyL2 = ihitL2 > 0 ? imatchL2 / ihitL2 : 0.0;
594 
595  m_monObj->setVariable("efficiency", var_efficiency);
596  m_monObj->setVariable("efficiencyL1", var_efficiencyL1);
597  m_monObj->setVariable("efficiencyL2", var_efficiencyL2);
598  m_monObj->setVariable("nmodules", ieff);
599 
600  setEpicsPV("Status", stat_data);
601  // only update if statistics is reasonable, we dont want "0" drops between runs!
602  if (stat_data != c_StatusTooFew) {
603  setEpicsPV("All", var_efficiency);
604  setEpicsPV("L1", var_efficiencyL1);
605  setEpicsPV("L2", var_efficiencyL2);
606  }
607  }
608 }
609 
611 {
612  B2DEBUG(1, "DQMHistAnalysisPXDEff: terminate called");
613 }
614 
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
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.
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, 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.