Belle II Software  release-05-01-25
DQMHistAnalysisSVDEfficiency.cc
1 //+
2 // File : DQMHistAnalysisSVDEfficiency.cc
3 // Description :
4 //
5 // Author : Giulia Casarosa (PI), Gaetano De Marino (PI)
6 // Date : 20190428
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisSVDEfficiency.h>
11 #include <vxd/geometry/GeoCache.h>
12 
13 #include <TROOT.h>
14 #include <TStyle.h>
15 #include <TString.h>
16 
17 using namespace std;
18 using namespace Belle2;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisSVDEfficiency)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
31 {
32  //Parameter definition
33  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: Constructor done.");
34 
35  addParam("RefHistoFile", m_refFileName, "Reference histrogram file name", std::string("SVDrefHisto.root"));
36  addParam("effLevel_Error", m_effError, "Efficiency error (%) level (red)", float(0.9));
37  addParam("effLevel_Warning", m_effWarning, "Efficiency WARNING (%) level (orange)", float(0.94));
38  addParam("effLevel_Empty", m_effEmpty, "Threshold to consider the sensor efficiency as too low", float(0));
39  addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(false));
40 }
41 
42 DQMHistAnalysisSVDEfficiencyModule::~DQMHistAnalysisSVDEfficiencyModule() { }
43 
44 void DQMHistAnalysisSVDEfficiencyModule::initialize()
45 {
46  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: initialized.");
47  B2DEBUG(10, " black = " << kBlack);
48  B2DEBUG(10, " green = " << kGreen);
49  B2DEBUG(10, " orange = " << kOrange);
50  B2DEBUG(10, " Red = " << kRed);
51 
52  m_refFile = NULL;
53  if (m_refFileName != "") {
54  m_refFile = new TFile(m_refFileName.data(), "READ");
55  }
56 
57  //search for reference
58  if (m_refFile && m_refFile->IsOpen()) {
59  B2INFO("SVD DQMHistAnalysis: reference root file (" << m_refFileName << ") FOUND, reading ref histograms");
60 
61  TH1F* ref_eff = (TH1F*)m_refFile->Get("refEfficiency");
62  if (!ref_eff)
63  B2WARNING("SVD DQMHistAnalysis: Efficiency Level Refence not found! using module parameters");
64  else {
65  m_effEmpty = ref_eff->GetBinContent(1);
66  m_effWarning = ref_eff->GetBinContent(2);
67  m_effError = ref_eff->GetBinContent(3);
68  }
69 
70  } else
71  B2WARNING("SVD DQMHistAnalysis: reference root file (" << m_refFileName << ") not found, or closed, using module parameters");
72 
73  B2INFO(" SVD efficiency thresholds:");
74  B2INFO(" EFFICIENCY: empty < " << m_effEmpty << " < normal < " << m_effWarning << " < warning < " << m_effError << " < error");
75 
76 
77  //build the legend
78  m_legProblem = new TPaveText(11, findBinY(4, 3) - 3, 16, findBinY(4, 3));
79  m_legProblem->AddText("ERROR!");
80  m_legProblem->AddText("at least one sensor with:");
81  m_legProblem->AddText(Form("efficiency < %1.0f%%", m_effError * 100));
82  m_legProblem->SetFillColor(kRed);
83  m_legWarning = new TPaveText(11, findBinY(4, 3) - 3, 16, findBinY(4, 3));
84  m_legWarning->AddText("WARNING!");
85  m_legWarning->AddText("at least one sensor with:");
86  m_legWarning->AddText(Form("%1.0f%% < efficiency < %1.0f%%", m_effError * 100, m_effWarning * 100));
87  m_legWarning->SetFillColor(kOrange);
88  m_legNormal = new TPaveText(11, findBinY(4, 3) - 3, 16, findBinY(4, 3));
89  m_legNormal->AddText("EFFICIENCY WITHIN LIMITS");
90  m_legNormal->AddText(Form("efficiency > %1.0f%%", m_effWarning * 100));
91  m_legNormal->SetFillColor(kGreen);
92  m_legNormal->SetBorderSize(0.);
93  m_legNormal->SetLineColor(kBlack);
94  m_legEmpty = new TPaveText(11, findBinY(4, 3) - 2, 16, findBinY(4, 3));
95  m_legEmpty->AddText("NO DATA RECEIVED");
96  m_legEmpty->AddText("from at least one sensor");
97  m_legEmpty->SetFillColor(kBlack);
98  m_legEmpty->SetTextColor(kWhite);
99  m_legEmpty->SetBorderSize(0.);
100  m_legEmpty->SetLineColor(kBlack);
101 
102 
103  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
104 
105  //collect the list of all SVD Modules in the geometry here
106  std::vector<VxdID> sensors = geo.getListOfSensors();
107  for (VxdID& aVxdID : sensors) {
108  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
109  // B2INFO("VXD " << aVxdID);
110  if (info.getType() != VXD::SensorInfoBase::SVD) continue;
111  m_SVDModules.push_back(aVxdID); // reorder, sort would be better
112  }
113  std::sort(m_SVDModules.begin(), m_SVDModules.end()); // back to natural order
114 
115  gROOT->cd();
116  m_cEfficiencyU = new TCanvas("SVDAnalysis/c_SVDEfficiencyU");
117  m_cEfficiencyV = new TCanvas("SVDAnalysis/c_SVDEfficiencyV");
118  m_cEfficiencyErrU = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrU");
119  m_cEfficiencyErrV = new TCanvas("SVDAnalysis/c_SVDEfficiencyErrV");
120 
121  m_hEfficiency = new SVDSummaryPlots("SVDEfficiency@view", "Summary of SVD efficiencies (%), @view/@side Side");
122  m_hEfficiencyErr = new SVDSummaryPlots("SVDEfficiencyErr@view", "Summary of SVD efficiencies errors (%), @view/@side Side");
123 }
124 
125 void DQMHistAnalysisSVDEfficiencyModule::beginRun()
126 {
127  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: beginRun called.");
128  m_cEfficiencyU->Clear();
129  m_cEfficiencyV->Clear();
130  m_cEfficiencyErrU->Clear();
131  m_cEfficiencyErrV->Clear();
132 }
133 
134 void DQMHistAnalysisSVDEfficiencyModule::event()
135 {
136  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: event called.");
137 
138  //SETUP gSTYLE - all plots
139  // gStyle->SetOptStat(0);
140  // gStyle->SetTitleY(.97);
141 
142  //check MODULE EFFICIENCY
143  m_effUstatus = 0;
144  m_effVstatus = 0;
145  m_effUErrstatus = 0;
146  m_effVErrstatus = 0;
147 
148  //set dedicate gStyle
149  // const Int_t colNum = 4;
150  // Int_t palette[colNum] {kBlack, kGreen, kOrange, kRed};
151  // gStyle->SetPalette(colNum, palette);
152  gStyle->SetOptStat(0);
153  gStyle->SetPaintTextFormat("2.1f");
154 
155  m_hEfficiency->getHistogram(0)->Reset();
156  m_hEfficiency->getHistogram(1)->Reset();
157  m_hEfficiencyErr->getHistogram(0)->Reset();
158  m_hEfficiencyErr->getHistogram(1)->Reset();
159 
160  Float_t effU;
161  Float_t effV;
162  Float_t erreffU;
163  Float_t erreffV;
164 
165  //Efficiency for the U side
166  TH2F* found_tracksU = (TH2F*)findHist("SVDEfficiency/TrackHitsU");
167  TH2F* matched_clusU = (TH2F*)findHist("SVDEfficiency/MatchedHitsU");
168 
169  if (matched_clusU == NULL || found_tracksU == NULL) {
170  B2INFO("Histograms needed for Efficiency computation are not found");
171  m_cEfficiencyU->SetFillColor(kRed);
172  } else {
173  B2INFO("U-side Before loop on sensors, size :" << m_SVDModules.size());
174  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
175  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
176  int bin = found_tracksU->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
177  m_SVDModules[i].getSensorNumber()));
178  float numU = matched_clusU->GetBinContent(bin);
179  float denU = found_tracksU->GetBinContent(bin);
180  if (denU > 0)
181  effU = numU / denU;
182  else
183  effU = -1;
184  B2DEBUG(10, "effU = " << numU << "/" << denU << " = " << effU);
185  m_hEfficiency->fill(m_SVDModules[i], 1, effU * 100);
186  if (effU == -1)
187  erreffU = -1;
188  else
189  erreffU = std::sqrt(effU * (1 - effU) / denU);
190  m_hEfficiencyErr->fill(m_SVDModules[i], 1, erreffU * 100);
191 
192  if (effU <= m_effEmpty) {
193  if (m_effUstatus < 1) m_effUstatus = 1;
194  } else if (effU < m_effWarning) {
195  if (effU > m_effError) {
196  if (m_effUstatus < 2) m_effUstatus = 2;
197  } else {
198  if (m_effUstatus < 3) m_effUstatus = 3;
199  }
200  }
201  }
202  }
203 
204  //Efficiency for the V side
205  TH2F* found_tracksV = (TH2F*)findHist("SVDEfficiency/TrackHitsV");
206  TH2F* matched_clusV = (TH2F*)findHist("SVDEfficiency/MatchedHitsV");
207 
208  if (matched_clusV == NULL || found_tracksV == NULL) {
209  B2INFO("Histograms needed for Efficiency computation are not found");
210  m_cEfficiencyV->SetFillColor(kRed);
211  } else {
212  B2INFO("V-side Before loop on sensors, size :" << m_SVDModules.size());
213  for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
214  B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
215  int bin = found_tracksV->FindBin(m_SVDModules[i].getLadderNumber(), findBinY(m_SVDModules[i].getLayerNumber(),
216  m_SVDModules[i].getSensorNumber()));
217  float numV = matched_clusV->GetBinContent(bin);
218  float denV = found_tracksV->GetBinContent(bin);
219  if (denV > 0)
220  effV = numV / denV;
221  else
222  effV = -1;
223 
224  B2DEBUG(10, "effV = " << numV << "/" << denV << " = " << effV);
225  m_hEfficiency->fill(m_SVDModules[i], 0, effV * 100);
226  if (effV == -1)
227  erreffV = -1;
228  else
229  erreffV = std::sqrt(effV * (1 - effV) / denV);
230 
231  m_hEfficiencyErr->fill(m_SVDModules[i], 0, erreffV * 100);
232 
233  if (effV <= m_effEmpty) {
234  if (m_effVstatus < 1) m_effVstatus = 1;
235  } else if (effV < m_effWarning) {
236  if (effV > m_effError) {
237  if (m_effVstatus < 2) m_effVstatus = 2;
238  } else {
239  if (m_effVstatus < 3) m_effVstatus = 3;
240  }
241  }
242  }
243  }
244 
245 
246  //update summary
247  m_cEfficiencyU->cd();
248 
249 
250 
251  m_hEfficiency->getHistogram(1)->Draw("text");
252 
253  if (m_effUstatus == 0) {
254  m_cEfficiencyU->SetFillColor(kGreen);
255  m_cEfficiencyU->SetFrameFillColor(10);
256  m_legNormal->Draw("same");
257  } else {
258  if (m_effUstatus == 3) {
259  m_cEfficiencyU->SetFillColor(kRed);
260  m_cEfficiencyU->SetFrameFillColor(10);
261  m_legProblem->Draw("same");
262  }
263  if (m_effUstatus == 2) {
264  m_cEfficiencyU->SetFillColor(kOrange);
265  m_cEfficiencyU->SetFrameFillColor(10);
266  m_legWarning->Draw("same");
267  }
268  if (m_effUstatus == 1) {
269  m_cEfficiencyU->SetFillColor(kGray);
270  m_cEfficiencyU->SetFrameFillColor(10);
271  m_legEmpty->Draw("same");
272  }
273  }
274  m_cEfficiencyU->Draw("text");
275  m_cEfficiencyU->Update();
276  m_cEfficiencyU->Modified();
277  m_cEfficiencyU->Update();
278 
279 
280  //update summary
281  m_cEfficiencyV->cd();
282  m_hEfficiency->getHistogram(0)->Draw("text");
283 
284  if (m_effVstatus == 0) {
285  m_cEfficiencyV->SetFillColor(kGreen);
286  m_cEfficiencyV->SetFrameFillColor(10);
287  m_legNormal->Draw("same");
288  } else {
289  if (m_effVstatus == 3) {
290  m_cEfficiencyV->SetFillColor(kRed);
291  m_cEfficiencyV->SetFrameFillColor(10);
292  m_legProblem->Draw("same");
293  }
294  if (m_effVstatus == 2) {
295  m_cEfficiencyV->SetFillColor(kOrange);
296  m_cEfficiencyV->SetFrameFillColor(10);
297  m_legWarning->Draw("same");
298  }
299  if (m_effVstatus == 1) {
300  m_cEfficiencyV->SetFillColor(kGray);
301  m_cEfficiencyV->SetFrameFillColor(10);
302  m_legEmpty->Draw("same");
303  }
304  }
305 
306  m_cEfficiencyV->Draw();
307  m_cEfficiencyV->Update();
308  m_cEfficiencyV->Modified();
309  m_cEfficiencyV->Update();
310 
311  m_cEfficiencyErrU->cd();
312  m_hEfficiencyErr->getHistogram(1)->Draw("colztext");
313  m_cEfficiencyErrU->Draw();
314  m_cEfficiencyErrU->Update();
315  m_cEfficiencyErrU->Modified();
316  m_cEfficiencyErrU->Update();
317 
318  m_cEfficiencyErrV->cd();
319  m_hEfficiencyErr->getHistogram(0)->Draw("colztext");
320  m_cEfficiencyErrV->Draw();
321  m_cEfficiencyErrV->Update();
322  m_cEfficiencyErrV->Modified();
323  m_cEfficiencyErrV->Update();
324 
325 
326 }
327 
328 void DQMHistAnalysisSVDEfficiencyModule::endRun()
329 {
330  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: endRun called");
331  if (m_printCanvas) {
332  m_cEfficiencyU->Print("c_SVDEfficiencyU.pdf");
333  m_cEfficiencyV->Print("c_SVDEfficiencyV.pdf");
334  m_cEfficiencyErrU->Print("c_SVDEfficiencyErrU.pdf");
335  m_cEfficiencyErrV->Print("c_SVDEfficiencyErrV.pdf");
336  }
337 }
338 
339 void DQMHistAnalysisSVDEfficiencyModule::terminate()
340 {
341  B2DEBUG(10, "DQMHistAnalysisSVDEfficiency: terminate called");
342 
343  delete m_refFile;
344  delete m_legProblem;
345  delete m_legWarning;
346  delete m_legNormal;
347  delete m_legEmpty;
348  delete m_hEfficiency;
349  delete m_cEfficiencyU;
350  delete m_cEfficiencyV;
351  delete m_hEfficiencyErr;
352  delete m_cEfficiencyErrU;
353  delete m_cEfficiencyErrV;
354 }
355 
356 // return y coordinate in TH2F histogram for specified sensor
357 Int_t DQMHistAnalysisSVDEfficiencyModule::findBinY(Int_t layer, Int_t sensor)
358 {
359  if (layer == 3)
360  return sensor; //2 -> 1,2
361  if (layer == 4)
362  return 2 + 1 + sensor; //6 -> 4,5,6
363  if (layer == 5)
364  return 6 + 1 + sensor; // 11 -> 8, 9, 10, 11
365  if (layer == 6)
366  return 11 + 1 + sensor; // 17 -> 13, 14, 15, 16, 17
367  else
368  return -1;
369 }
370 
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::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
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisSVDEfficiencyModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisSVDEfficiency.h:25
Belle2::SVDSummaryPlots
class to summarize SVD quantities per sensor and side
Definition: SVDSummaryPlots.h:35
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::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27