Belle II Software  release-08-01-10
DQMHistAnalysisSVDEfficiency.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 : DQMHistAnalysisSVDEfficiency.cc
10 // Description : module for DQM histogram analysis of SVD sensors efficiencies
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisSVDEfficiency.h>
15 #include <vxd/geometry/GeoCache.h>
16 
17 #include <TROOT.h>
18 #include <TStyle.h>
19 #include <TString.h>
20 
21 using namespace std;
22 using namespace Belle2;
23 
24 //-----------------------------------------------------------------
25 // Register the Module
26 //-----------------------------------------------------------------
27 REG_MODULE(DQMHistAnalysisSVDEfficiency);
28 
29 //-----------------------------------------------------------------
30 // Implementation
31 //-----------------------------------------------------------------
32 
33 DQMHistAnalysisSVDEfficiencyModule::DQMHistAnalysisSVDEfficiencyModule()
35  m_effUstatus(good),
36  m_effVstatus(good)
37 {
38  //Parameter definition
39  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: Constructor done.");
40 
41  setDescription("DQM Analysis Module that computes the average SVD sensor efficiency.");
42 
43  addParam("RefHistoFile", m_refFileName, "Reference histogram file name", std::string("SVDrefHisto.root"));
44  addParam("effLevel_Error", m_effError, "Efficiency error (%) level (red)", double(0.9));
45  addParam("effLevel_Warning", m_effWarning, "Efficiency WARNING (%) level (orange)", double(0.94));
46  addParam("statThreshold", m_statThreshold, "minimal number of tracks per sensor to set green/red alert", double(100));
47  addParam("samples3", m_3Samples, "if True 3 samples histograms analysis is performed", bool(false));
48  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("SVD:"));
49 }
50 
52 
54 {
55  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: initialize");
56 
57  //build the legend
58  m_legProblem = new TPaveText(0.62, 0.22, 0.88, 0.35, "brNDC");
59  m_legProblem->AddText("ERROR!");
60  m_legProblem->AddText("at least one sensor with:");
61  m_legProblem->AddText(Form("efficiency < %1.0f%%", m_effError * 100));
62  m_legProblem->SetFillColor(c_ColorDefault);
63  m_legProblem->SetLineColor(kBlack);
64 
65  m_legWarning = new TPaveText(0.62, 0.22, 0.88, 0.35, "brNDC");
66  m_legWarning->AddText("WARNING!");
67  m_legWarning->AddText("at least one sensor with:");
68  m_legWarning->AddText(Form("%1.0f%% < efficiency < %1.0f%%", m_effError * 100, m_effWarning * 100));
69  m_legWarning->SetFillColor(c_ColorDefault);
70  m_legWarning->SetLineColor(kBlack);
71 
72  m_legNormal = new TPaveText(0.62, 0.22, 0.88, 0.35, "brNDC");
73  m_legNormal->AddText("EFFICIENCY WITHIN LIMITS");
74  m_legNormal->AddText(Form("efficiency > %1.0f%%", m_effWarning * 100));
75  m_legNormal->SetFillColor(c_ColorDefault);
76  m_legNormal->SetLineColor(kBlack);
77 
78  m_legEmpty = new TPaveText(0.62, 0.22, 0.88, 0.35, "brNDC");
79  m_legEmpty->AddText("Not enough statistics,");
80  m_legEmpty->AddText("check again in a few minutes");
81  m_legEmpty->SetFillColor(c_ColorDefault);
82  m_legEmpty->SetTextColor(kWhite);
83  m_legEmpty->SetLineColor(kBlack);
84 
85 
87 
88  //collect the list of all SVD Modules in the geometry here
89  std::vector<VxdID> sensors = geo.getListOfSensors();
90  for (VxdID& aVxdID : sensors) {
91  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
92  // B2INFO("VXD " << aVxdID);
93  if (info.getType() != VXD::SensorInfoBase::SVD) continue;
94  m_SVDModules.push_back(aVxdID); // reorder, sort would be better
95  }
96  std::sort(m_SVDModules.begin(), m_SVDModules.end()); // back to natural order
97 
98 
99  gROOT->cd();
100  m_cEfficiencyU = new TCanvas("SVDAnalysis/c_SVDEfficiencyU");
101  m_cEfficiencyV = new TCanvas("SVDAnalysis/c_SVDEfficiencyV");
102  m_cEfficiencyErrU = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrU");
103  m_cEfficiencyErrV = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrV");
104 
105  m_hEfficiency = new SVDSummaryPlots("SVDEfficiency@view", "Summary of SVD efficiencies (%), @view/@side Side");
107  m_hEfficiencyErr = new SVDSummaryPlots("SVDEfficiencyErr@view", "Summary of SVD efficiencies errors (%), @view/@side Side");
109 
110 
111  if (m_3Samples) {
112  m_cEfficiencyU3Samples = new TCanvas("SVDAnalysis/c_SVDEfficiencyU3Samples");
113  m_cEfficiencyV3Samples = new TCanvas("SVDAnalysis/c_SVDEfficiencyV3Samples");
114  m_cEfficiencyErrU3Samples = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrU3Samples");
115  m_cEfficiencyErrV3Samples = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrV3Samples");
116 
117  m_hEfficiency3Samples = new SVDSummaryPlots("SVD3Efficiency@view",
118  "Summary of SVD efficiencies (%), @view/@side Side for 3 samples");
120 
121  m_hEfficiencyErr3Samples = new SVDSummaryPlots("SVD3EfficiencyErr@view",
122  "Summary of SVD efficiencies errors (%), @view/@side Side for 3 samples");
124 
125  }
126 
127  //register limits for EPICS
128  registerEpicsPV(m_pvPrefix + "efficiencyLimits", "effLimits");
129 
130  //find nEvents testing if histograms are present
131  TH1* hnEvnts = findHist("SVDExpReco/SVDDQM_nEvents");
132  if (hnEvnts == NULL) {
133  B2INFO("no events, nothing to do here");
134  return;
135  }
136 
137 }
138 
140 {
141  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: beginRun called.");
142 
143  if (m_cEfficiencyU)
144  m_cEfficiencyU->Clear();
145  if (m_cEfficiencyV)
146  m_cEfficiencyV->Clear();
147  if (m_cEfficiencyErrU)
148  m_cEfficiencyErrU->Clear();
149  if (m_cEfficiencyErrV)
150  m_cEfficiencyErrV->Clear();
151 
152  if (m_3Samples) {
154  m_cEfficiencyU3Samples->Clear();
156  m_cEfficiencyV3Samples->Clear();
158  m_cEfficiencyErrU3Samples->Clear();
160  m_cEfficiencyErrV3Samples->Clear();
161  }
162 
163  //Retrieve limits from EPICS
164  double effErrorLo = 0.;
165  double effWarnLo = 0.;
166 
167  requestLimitsFromEpicsPVs("effLimits", effErrorLo, effWarnLo, m_effWarning, m_effError);
168 
169  B2DEBUG(10, " SVD efficiency thresholds taken from EPICS configuration file:");
170  B2DEBUG(10, " EFFICIENCY: normal > " << m_effWarning << " > warning > " << m_effError << " > error with minimum statistics of " <<
172 }
173 
175 {
176  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: event called.");
177 
178  //find nEvents
179  TH1* hnEvnts = findHist("SVDExpReco/SVDDQM_nEvents", true);
180  if (hnEvnts == NULL) {
181  B2INFO("no events, nothing to do here");
182  return;
183  } else {
184  B2DEBUG(10, "SVDExpReco/SVDDQM_nEvents found");
185  }
186 
187  TString tmp = hnEvnts->GetTitle();
188  Int_t pos = tmp.Last('~');
189  if (pos == -1) pos = 0;
190 
191  TString runID = tmp(pos, tmp.Length() - pos);
192  B2INFO("DQMHistAnalysisSVDEfficiencyModule::runID = " << runID);
193 
194  gStyle->SetOptStat(0);
195  gStyle->SetPaintTextFormat("2.1f");
196 
197  // do it by nhand, the interface of the SVDSummaryPlots does not allow to change the title after cstr
198  if (m_hEfficiency) {
199  m_hEfficiency->reset();
200  m_hEfficiency->setRunID(runID);
201  }
202 
203  if (m_hEfficiencyErr) {
205  m_hEfficiencyErr->setRunID(runID);
206  }
207 
208  if (m_3Samples) {
209  if (m_hEfficiency3Samples) {
212  }
213 
217  }
218  }
219 
220  Float_t effU = -1;
221  Float_t effV = -1;
222  Float_t erreffU = -1;
223  Float_t erreffV = -1;
224 
225  // Efficiency for the U side
226  TH2F* found_tracksU = (TH2F*)findHist("SVDEfficiency/TrackHitsU");
227  TH2F* matched_clusU = (TH2F*)findHist("SVDEfficiency/MatchedHitsU");
228 
229  if (matched_clusU == NULL || found_tracksU == NULL) {
230  B2INFO("Histograms needed for Efficiency computation are not found");
231  m_cEfficiencyU->Draw();
232  m_cEfficiencyU->cd();
233  m_hEfficiency->getHistogram(1)->Draw("text");
235  } else {
236  B2DEBUG(10, "U-side Before loop on sensors, size :" << m_SVDModules.size());
237  m_effUstatus = good;
238  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
239  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
240  int bin = found_tracksU->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
241  m_SVDModules[i].getSensorNumber()));
242  float numU = matched_clusU->GetBinContent(bin);
243  float denU = found_tracksU->GetBinContent(bin);
244  if (denU > 0)
245  effU = numU / denU;
246  else
247  effU = -1;
248  B2DEBUG(10, "effU = " << numU << "/" << denU << " = " << effU);
249  m_hEfficiency->fill(m_SVDModules[i], 1, effU * 100);
250  if (effU == -1)
251  erreffU = -1;
252  else
253  erreffU = std::sqrt(effU * (1 - effU) / denU);
254  m_hEfficiencyErr->fill(m_SVDModules[i], 1, erreffU * 100);
255 
256  if (denU < m_statThreshold) {
257  m_effUstatus = std::max(lowStat, m_effUstatus);
258  } else if (effU > m_effWarning) {
259  m_effUstatus = std::max(good, m_effUstatus);
260  } else if ((effU <= m_effWarning) && (effU > m_effError)) {
261  m_effUstatus = std::max(warning, m_effUstatus);
262  } else if ((effU <= m_effError)) {
263  m_effUstatus = std::max(error, m_effUstatus);
264  }
265  B2DEBUG(10, "Status is " << m_effUstatus);
266  }
267  }
268 
269  //Efficiency for the V side
270  TH2F* found_tracksV = (TH2F*)findHist("SVDEfficiency/TrackHitsV");
271  TH2F* matched_clusV = (TH2F*)findHist("SVDEfficiency/MatchedHitsV");
272 
273  if (matched_clusV == NULL || found_tracksV == NULL) {
274  B2INFO("Histograms needed for Efficiency computation are not found");
275  m_cEfficiencyV->cd();
276  m_cEfficiencyV->Draw();
277  m_hEfficiency->getHistogram(0)->Draw("text");
279  } else {
280  B2DEBUG(10, "V-side Before loop on sensors, size :" << m_SVDModules.size());
281  m_effVstatus = good;
282  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
283  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
284  int bin = found_tracksV->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
285  m_SVDModules[i].getSensorNumber()));
286  float numV = matched_clusV->GetBinContent(bin);
287  float denV = found_tracksV->GetBinContent(bin);
288  if (denV > 0)
289  effV = numV / denV;
290  else
291  effV = -1;
292 
293  B2DEBUG(10, "effV = " << numV << "/" << denV << " = " << effV);
294  m_hEfficiency->fill(m_SVDModules[i], 0, effV * 100);
295  if (effV == -1)
296  erreffV = -1;
297  else
298  erreffV = std::sqrt(effV * (1 - effV) / denV);
299 
300  m_hEfficiencyErr->fill(m_SVDModules[i], 0, erreffV * 100);
301 
302  if (denV < m_statThreshold) {
303  m_effVstatus = std::max(lowStat, m_effVstatus);
304  } else if (effV > m_effWarning) {
305  m_effVstatus = std::max(good, m_effVstatus);
306  } else if ((effV <= m_effWarning) && (effV > m_effError)) {
307  m_effVstatus = std::max(warning, m_effVstatus);
308  } else if ((effV <= m_effError)) {
309  m_effVstatus = std::max(error, m_effVstatus);
310  }
311  B2DEBUG(10, "Status is " << m_effVstatus);
312  }
313  }
314 
315  // update summary for U side
316  m_cEfficiencyU->Draw();
317  m_cEfficiencyU->cd();
318  if (m_hEfficiency)
319  m_hEfficiency->getHistogram(1)->Draw("text");
320 
321  switch (m_effUstatus) {
322  case good: {
324  m_legNormal->Draw();
325  break;
326  }
327  case error: {
329  m_legProblem->Draw();
330  break;
331  }
332  case warning: {
334  m_legWarning->Draw();
335  break;
336  }
337  case lowStat: {
339  m_legEmpty->Draw();
340  break;
341  }
342  default: {
343  B2INFO("effUstatus not set properly: " << m_effUstatus);
344  break;
345  }
346  }
347 // setEpicsPV("EfficiencyUAlarm", alarm);
348 
349  m_cEfficiencyU->Update();
350  m_cEfficiencyU->Modified();
351  m_cEfficiencyU->Update();
352 
353  // update summary for V side
354  m_cEfficiencyV->cd();
355  m_cEfficiencyV->Draw();
356  if (m_hEfficiency)
357  m_hEfficiency->getHistogram(0)->Draw("text");
358 
359  switch (m_effVstatus) {
360  case good: {
362  m_legNormal->Draw();
363  break;
364  }
365  case error: {
367  m_legProblem->Draw();
368  break;
369  }
370  case warning: {
372  m_legWarning->Draw();
373  break;
374  }
375  case lowStat: {
377  m_legEmpty->Draw();
378  break;
379  }
380  default: {
381  B2INFO("effVstatus not set properly: " << m_effVstatus);
382  break;
383  }
384  }
385 
386  m_cEfficiencyV->Update();
387  m_cEfficiencyV->Modified();
388  m_cEfficiencyV->Update();
389 
390  m_cEfficiencyErrU->cd();
391  if (m_hEfficiencyErr)
392  m_hEfficiencyErr->getHistogram(1)->Draw("colztext");
393  m_cEfficiencyErrU->Draw();
394  m_cEfficiencyErrU->Update();
395  m_cEfficiencyErrU->Modified();
396  m_cEfficiencyErrU->Update();
397 
398  m_cEfficiencyErrV->cd();
399  if (m_hEfficiencyErr)
400  m_hEfficiencyErr->getHistogram(0)->Draw("colztext");
401  m_cEfficiencyErrV->Draw();
402  m_cEfficiencyErrV->Update();
403  m_cEfficiencyErrV->Modified();
404  m_cEfficiencyErrV->Update();
405 
406  if (m_3Samples) {
408  m_hEfficiency3Samples->getHistogram(0)->Reset();
409  m_hEfficiency3Samples->getHistogram(1)->Reset();
412 
413  // Efficiency for the U side - 3 samples
414  TH2F* found3_tracksU = (TH2F*)findHist("SVDEfficiency/TrackHits3U");
415  TH2F* matched3_clusU = (TH2F*)findHist("SVDEfficiency/MatchedHits3U");
416 
417  if (matched3_clusU == NULL || found3_tracksU == NULL) {
418  B2INFO("Histograms needed for Efficiency computation are not found");
419  m_cEfficiencyU3Samples->Draw();
421  m_hEfficiency3Samples->getHistogram(1)->Draw("text");
423  } else {
424  B2DEBUG(10, "U-side Before loop on sensors, size :" << m_SVDModules.size());
425  m_effUstatus = good;
426  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
427  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
428  int bin = found3_tracksU->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
429  m_SVDModules[i].getSensorNumber()));
430  float numU = matched3_clusU->GetBinContent(bin);
431  float denU = found3_tracksU->GetBinContent(bin);
432  if (denU > 0)
433  effU = numU / denU;
434  else
435  effU = -1;
436  B2DEBUG(10, "effU = " << numU << "/" << denU << " = " << effU);
437 
438  m_hEfficiency3Samples->fill(m_SVDModules[i], 1, effU * 100);
439  if (effU == -1)
440  erreffU = -1;
441  else
442  erreffU = std::sqrt(effU * (1 - effU) / denU);
443  m_hEfficiencyErr3Samples->fill(m_SVDModules[i], 1, erreffU * 100);
444 
445  if (denU < m_statThreshold) {
446  m_effUstatus = std::max(lowStat, m_effUstatus);
447  } else if (effU > m_effWarning) {
448  m_effUstatus = std::max(good, m_effUstatus);
449  } else if ((effU <= m_effWarning) && (effU > m_effError)) {
450  m_effUstatus = std::max(warning, m_effUstatus);
451  } else if ((effU <= m_effError)) {
452  m_effUstatus = std::max(error, m_effUstatus);
453  }
454  B2DEBUG(10, "Status is " << m_effUstatus);
455  }
456  }
457 
458  //Efficiency for the V side - 3 samples
459  TH2F* found3_tracksV = (TH2F*)findHist("SVDEfficiency/TrackHits3V");
460  TH2F* matched3_clusV = (TH2F*)findHist("SVDEfficiency/MatchedHits3V");
461 
462  if (matched3_clusV == NULL || found3_tracksV == NULL) {
463  B2INFO("Histograms needed for Efficiency computation are not found");
464  m_cEfficiencyV3Samples->Draw();
466  m_hEfficiency3Samples->getHistogram(1)->Draw("text");
468  } else {
469  B2DEBUG(10, "V-side Before loop on sensors, size :" << m_SVDModules.size());
470  m_effVstatus = good;
471  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
472  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
473  int bin = found3_tracksV->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
474  m_SVDModules[i].getSensorNumber()));
475  float numV = matched3_clusV->GetBinContent(bin);
476  float denV = found3_tracksV->GetBinContent(bin);
477  if (denV > 0)
478  effV = numV / denV;
479  else
480  effV = -1;
481 
482  B2DEBUG(10, "effV = " << numV << "/" << denV << " = " << effV);
483  m_hEfficiency3Samples->fill(m_SVDModules[i], 0, effV * 100);
484  if (effV == -1)
485  erreffV = -1;
486  else
487  erreffV = std::sqrt(effV * (1 - effV) / denV);
488 
489  m_hEfficiencyErr3Samples->fill(m_SVDModules[i], 0, erreffV * 100);
490 
491  if (denV < m_statThreshold) {
492  m_effVstatus = std::max(lowStat, m_effVstatus);
493  } else if (effV > m_effWarning) {
494  m_effVstatus = std::max(good, m_effVstatus);
495  } else if ((effV <= m_effWarning) && (effV > m_effError)) {
496  m_effVstatus = std::max(warning, m_effVstatus);
497  } else if ((effV <= m_effError)) {
498  m_effVstatus = std::max(error, m_effVstatus);
499  }
500  B2DEBUG(10, "Status is " << m_effVstatus);
501  }
502  }
503 
504  // update summary for U side
505  m_cEfficiencyU3Samples->Draw();
508  m_hEfficiency3Samples->getHistogram(1)->Draw("text");
509 
510  switch (m_effUstatus) {
511  case good: {
513  m_legNormal->Draw();
514  break;
515  }
516  case error: {
518  m_legProblem->Draw();
519  break;
520  }
521  case warning: {
523  m_legWarning->Draw();
524  break;
525  }
526  case lowStat: {
528  m_legEmpty->Draw();
529  break;
530  }
531  default: {
532  B2INFO("effUstatus not set properly: " << m_effUstatus);
533  break;
534  }
535  }
536 
537  m_cEfficiencyU3Samples->Update();
538  m_cEfficiencyU3Samples->Modified();
539  m_cEfficiencyU3Samples->Update();
540 
541  // update summary for V side
542  m_cEfficiencyV3Samples->Draw();
545  m_hEfficiency3Samples->getHistogram(0)->Draw("text");
546 
547  switch (m_effVstatus) {
548  case good: {
550  m_legNormal->Draw();
551  break;
552  }
553  case error: {
555  m_legProblem->Draw();
556  break;
557  }
558  case warning: {
560  m_legWarning->Draw();
561  break;
562  }
563  case lowStat: {
565  m_legEmpty->Draw();
566  break;
567  }
568  default: {
569  B2INFO("effVstatus not set properly: " << m_effVstatus);
570  break;
571  }
572  }
573 
574  m_cEfficiencyV3Samples->Update();
575  m_cEfficiencyV3Samples->Modified();
576  m_cEfficiencyV3Samples->Update();
577 
580  m_hEfficiencyErr3Samples->getHistogram(1)->Draw("colztext");
582  m_cEfficiencyErrU3Samples->Update();
583  m_cEfficiencyErrU3Samples->Modified();
584  m_cEfficiencyErrU3Samples->Update();
585 
588  m_hEfficiencyErr3Samples->getHistogram(0)->Draw("colztext");
590  m_cEfficiencyErrV3Samples->Update();
591  m_cEfficiencyErrV3Samples->Modified();
592  m_cEfficiencyErrV3Samples->Update();
593  }
594 }
595 
597 {
598  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: endRun called");
599 }
600 
602 {
603  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: terminate called");
604 
605  delete m_refFile;
606  delete m_legProblem;
607  delete m_legWarning;
608  delete m_legNormal;
609  delete m_legEmpty;
610 
611  delete m_hEfficiency;
612  delete m_cEfficiencyU;
613  delete m_cEfficiencyV;
614  delete m_hEfficiencyErr;
615  delete m_cEfficiencyErrU;
616  delete m_cEfficiencyErrV;
617 
618  if (m_3Samples) {
619  delete m_hEfficiency3Samples;
620  delete m_cEfficiencyU3Samples;
621  delete m_cEfficiencyV3Samples;
625  }
626 }
627 
628 // return y coordinate in TH2F histogram for specified sensor
629 Int_t DQMHistAnalysisSVDEfficiencyModule::findBinY(Int_t layer, Int_t sensor)
630 {
631  if (layer == 3)
632  return sensor; //2 -> 1,2
633  if (layer == 4)
634  return 2 + 1 + sensor; //6 -> 4,5,6
635  if (layer == 5)
636  return 6 + 1 + sensor; // 11 -> 8, 9, 10, 11
637  if (layer == 6)
638  return 11 + 1 + sensor; // 17 -> 13, 14, 15, 16, 17
639  else
640  return -1;
641 }
642 
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).
@ c_ColorDefault
default for non-coloring
@ c_StatusDefault
default for non-coloring
@ c_StatusTooFew
Not enough entries/event to judge.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusWarning
Analysis result: Warning, there may be minor issues.
@ c_StatusGood
Analysis result: Good.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
SVDSummaryPlots * m_hEfficiency3Samples
efficiency histo for 3 samples
TCanvas * m_cEfficiencyErrV3Samples
efficiency V error plot canvas for 3 samples
TPaveText * m_legEmpty
efficiency plot legend, empty
double m_statThreshold
minimal number of tracks per sensor to set green or red frame
TCanvas * m_cEfficiencyErrU
efficiency U error plot canvas
Int_t findBinY(Int_t layer, Int_t sensor)
find Y bin corresponding to sensor, efficiency plot
std::vector< VxdID > m_SVDModules
IDs of all SVD Modules to iterate over.
effStatus m_effUstatus
number representing the status of the efficiency U side
SVDSummaryPlots * m_hEfficiencyErr3Samples
efficiency error histo for 3 samples
std::string m_pvPrefix
string prefix for EPICS PVs
double m_effWarning
warning level of the efficiency
double m_effError
error level of the efficiency
void terminate() override final
This method is called at the end of the event processing.
TPaveText * m_legWarning
efficiency plot legend, warning
void event() override final
This method is called for each event.
bool m_3Samples
if true enable 3 samples histograms analysis
TCanvas * m_cEfficiencyU3Samples
efficiency U plot canvas for 3 samples
void endRun() override final
This method is called if the current run ends.
TCanvas * m_cEfficiencyErrU3Samples
efficiency U error plot canvas for 3 samples
effStatus m_effVstatus
number representing the status of the efficiency V side
void beginRun() override final
Called when entering a new run.
TPaveText * m_legNormal
efficiency plot legend, normal
TCanvas * m_cEfficiencyV3Samples
efficiency V plot canvas for 3 samples
TCanvas * m_cEfficiencyErrV
efficiency V error plot canvas
TPaveText * m_legProblem
efficiency plot legend, problem
SVDSummaryPlots * m_hEfficiencyErr
efficiency error histo
TFile * m_refFile
The pointer to the reference file.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
class to summarize SVD quantities per sensor and side
void setStats(bool stats=true)
set histograms stat
void fill(int layer, int ladder, int sensor, int view, float value)
fill the histogram for
TH2F * getHistogram(int view)
get a reference to the histogram for
void setRunID(const TString &runID)
set run ids in title
void reset()
Reset histograms.
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.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.