Belle II Software  release-05-01-25
DQMHistAnalysisSVDOnMiraBelle.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2020 - Belle II Collaboration *
4  * *
5  * SVD Data Quality Monitor creates monitoring variables and pack them to *
6  * the Monitoring Objects, which are used as input to the Mirabelle *
7  * on the online system to show run by run dynamic of the variables. *
8  * This module provides higher level analyses on collected SVD DQM *
9  * histograms: *
10  * - SVDExpReco/SVDDQM_StripCountsU, SVDExpReco/SVDDQM_StripCountsV *
11  * - SVDExpReco/SVDDQM_nEvents *
12  * - SVDEfficiency/TrackHitsU, SVDEfficiency/MatchedHitsU *
13  * - SVDEfficiency/TrackHitsV, SVDEfficiency/MatchedHitsV *
14  * - SVDClsTrk/SVDTRK_ClusterChargeU3, SVDClsTrk/SVDTRK_ClusterChargeV3 *
15  * - SVDClsTrk/SVDTRK_ClusterChargeU456 *
16  * - SVDClsTrk/SVDTRK_ClusterChargeV456 *
17  * - SVDClsTrk/SVDTRK_ClusterSNRU3, SVDClsTrk/SVDTRK_ClusterSNRV3 *
18  * - SVDClsTrk/SVDTRK_ClusterSNRU456, SVDClsTrk/SVDTRK_ClusterSNRV456 *
19  * - SVDClsTrk/SVDTRK_ClusterTimeU3, SVDClsTrk/SVDTRK_ClusterTimeV3 *
20  * - SVDClsTrk/SVDTRK_ClusterSNRU456, SVDClsTrk/SVDTRK_ClusterSNRV456 *
21  * - SVDClsTrk/SVDTRK_StripMaxBinUAll, SVDClsTrk/SVDTRK_StripMaxBinVAll *
22  * *
23  * Author: The Belle II Collaboration *
24  * Contributors: Karol Adamczyk (Karol.Adamczyk@ifj.edu.pl) *
25  * *
26  * This software is provided "as is" without any warranty. *
27  **************************************************************************/
28 
29 #include <dqm/analysis/modules/DQMHistAnalysisSVDOnMiraBelle.h>
30 
31 using namespace std;
32 using namespace Belle2;
33 
34 //-----------------------------------------------------------------
35 // Register the Module
36 //-----------------------------------------------------------------
37 REG_MODULE(DQMHistAnalysisSVDOnMiraBelle)
38 
39 //-----------------------------------------------------------------
40 // Implementation
41 //-----------------------------------------------------------------
42 
45 {
46  setDescription("Extract monitoring variables from SVD DQM histograms");
47  setPropertyFlags(c_ParallelProcessingCertified);
48  B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: Constructor done.");
49 }
50 
51 DQMHistAnalysisSVDOnMiraBelleModule::~DQMHistAnalysisSVDOnMiraBelleModule() { }
52 
53 void DQMHistAnalysisSVDOnMiraBelleModule::initialize()
54 {
55  gROOT->cd();
56 
57  // add MonitoringObject
58  m_monObj = getMonitoringObject("svd");
59 
60  // list of canvases to be added to MonitoringObject
61  m_c_avgEfficiency = new TCanvas("svd_avgEfficiency", "matched clusters and found tracks", 0, 0, 800, 600);
62  m_c_avgOffOccupancy = new TCanvas("svd_avgOffOccupancy", "strips", 0, 0, 800, 600);
63  m_c_MPVChargeClusterOnTrack = new TCanvas("svd_MPVChargeClusterOnTrack", "charge from Clusters on Track Charge", 0, 0, 400, 400);
64  m_c_MPVSNRClusterOnTrack = new TCanvas("svd_MPVSNRClusterOnTrack", "SNR from Clusters on Track Charge", 0, 0, 400, 400);
65  m_c_MPVTimeClusterOnTrack = new TCanvas("svd_MPVTimeClusterOnTrack", "time from Clusters on Track Charge", 0, 0, 400, 400);
66  m_c_avgMaxBinClusterOnTrack = new TCanvas("svd_avgMaxBin", "average MaxBin", 0, 0, 800, 600);
67 
68  // add canvases used to create monitoring variables to MonitoringObject
69  m_monObj->addCanvas(m_c_avgEfficiency);
70  m_monObj->addCanvas(m_c_avgOffOccupancy);
71  m_monObj->addCanvas(m_c_MPVChargeClusterOnTrack);
72  m_monObj->addCanvas(m_c_MPVSNRClusterOnTrack);
73  m_monObj->addCanvas(m_c_MPVTimeClusterOnTrack);
74  m_monObj->addCanvas(m_c_avgMaxBinClusterOnTrack);
75 
76  B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: initialized.");
77 }
78 
79 void DQMHistAnalysisSVDOnMiraBelleModule::beginRun()
80 {
81  B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: beginRun called.");
82 }
83 
84 void DQMHistAnalysisSVDOnMiraBelleModule::event()
85 {
86  B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: event called.");
87 }
88 
89 void DQMHistAnalysisSVDOnMiraBelleModule::endRun()
90 {
91  // offline occupancy - integrated number of ZS5 fired strips
92  TH1F* h_zs5countsU = (TH1F*)findHist("SVDExpReco/SVDDQM_StripCountsU"); // made by SVDDQMExperssRecoModule
93  TH1F* h_zs5countsV = (TH1F*)findHist("SVDExpReco/SVDDQM_StripCountsV");
94  TH1F* h_events = (TH1F*)findHist("SVDExpReco/SVDDQM_nEvents");
95 
96  // adding histograms to canvas
97  m_c_avgOffOccupancy->Clear();
98  m_c_avgOffOccupancy->Divide(2, 2);
99  m_c_avgOffOccupancy->cd(1);
100  if (h_zs5countsU) h_zs5countsU->Draw("colz");
101  m_c_avgOffOccupancy->cd(2);
102  if (h_zs5countsV) h_zs5countsU->Draw("colz");
103  m_c_avgOffOccupancy->cd(3);
104  if (h_events) h_events->Draw("colz");
105 
106  int nE = h_events->GetEntries(); // number of events for all "clusters on track" histograms
107 
108  // average occupancy for each layer
109  std::vector<float> avgOffOccL3UV = avgOccupancyUV(3, h_zs5countsU, h_zs5countsV, 0, 14, 1, 1, nE);
110 
111  std::vector<float> avgOffOccL4UV = avgOccupancyUV(4, h_zs5countsU, h_zs5countsV, 0, 30, 15, 1, nE);
112 
113  std::vector<float> avgOffOccL5UV = avgOccupancyUV(5, h_zs5countsU, h_zs5countsV, 0, 48, 45, 1, nE);
114 
115  std::vector<float> avgOffOccL6UV = avgOccupancyUV(6, h_zs5countsU, h_zs5countsV, 0, 80, 93, 1, nE);
116 
117  // average occupancy for middle plane sensors
118  std::vector<float> avgOffOccL4X2UV = avgOccupancyUV(4, h_zs5countsU, h_zs5countsV, 0, 10, 16, 3, nE); // L4.X.2
119 
120  std::vector<float> avgOffOccL5X2UV = avgOccupancyUV(5, h_zs5countsU, h_zs5countsV, 0, 12, 36, 4, nE); // L5.X.2
121 
122  std::vector<float> avgOffOccL5X3UV = avgOccupancyUV(5, h_zs5countsU, h_zs5countsV, 0, 12, 37, 4, nE); // L5.X.3
123 
124  std::vector<float> avgOffOccL6X2UV = avgOccupancyUV(6, h_zs5countsU, h_zs5countsV, 0, 16, 94, 5, nE); // L6.X.2
125 
126  std::vector<float> avgOffOccL6X3UV = avgOccupancyUV(6, h_zs5countsU, h_zs5countsV, 0, 16, 95, 5, nE); // L6.X.3
127 
128  std::vector<float> avgOffOccL6X4UV = avgOccupancyUV(6, h_zs5countsU, h_zs5countsV, 0, 16, 96, 5, nE); // L6.X.4
129 
130  // average occupancy for high occupancy sensors
131  std::vector<float> avgOffOccL311UV = highOccupancySensor(3, h_zs5countsU, h_zs5countsV, 1, nE); // L3.1.1
132 
133  std::vector<float> avgOffOccL312UV = highOccupancySensor(3, h_zs5countsU, h_zs5countsV, 2, nE); // L3.1.2
134 
135  std::vector<float> avgOffOccL321UV = highOccupancySensor(3, h_zs5countsU, h_zs5countsV, 3, nE); // L3.2.1
136 
137  std::vector<float> avgOffOccL322UV = highOccupancySensor(3, h_zs5countsU, h_zs5countsV, 4, nE); // L3.2.2
138 
139  std::vector<float> avgOffOccL461UV = highOccupancySensor(4, h_zs5countsU, h_zs5countsV, 30, nE); // L4.6.1
140 
141  std::vector<float> avgOffOccL462UV = highOccupancySensor(4, h_zs5countsU, h_zs5countsV, 31, nE); // L4.6.2
142 
143  std::vector<float> avgOffOccL581UV = highOccupancySensor(5, h_zs5countsU, h_zs5countsV, 73, nE); // L5.8.1
144 
145  std::vector<float> avgOffOccL582UV = highOccupancySensor(5, h_zs5countsU, h_zs5countsV, 74, nE); // L5.8.2
146 
147  std::vector<float> avgOffOccL6101UV = highOccupancySensor(6, h_zs5countsU, h_zs5countsV, 138, nE); // L6.10.1
148 
149  std::vector<float> avgOffOccL6102UV = highOccupancySensor(6, h_zs5countsU, h_zs5countsV, 139, nE); // L6.10.2
150 
151  // setting monitoring variables
152  if (h_zs5countsU == NULL || h_events == NULL) {
153  B2INFO("Histograms needed for Average Offline Occupancy on U side are not found");
154  m_monObj->setVariable("avgOffOccL3U", -1);
155  m_monObj->setVariable("avgOffOccL4U", -1);
156  m_monObj->setVariable("avgOffOccL5U", -1);
157  m_monObj->setVariable("avgOffOccL6U", -1);
158  m_monObj->setVariable("avgOffOccL4X2U", -1);
159  m_monObj->setVariable("avgOffOccL5X2U", -1);
160  m_monObj->setVariable("avgOffOccL5X3U", -1);
161  m_monObj->setVariable("avgOffOccL6X2U", -1);
162  m_monObj->setVariable("avgOffOccL6X3U", -1);
163  m_monObj->setVariable("avgOffOccL6X4U", -1);
164  m_monObj->setVariable("avgOffOccL311U", -1);
165  m_monObj->setVariable("avgOffOccL312U", -1);
166  m_monObj->setVariable("avgOffOccL321U", -1);
167  m_monObj->setVariable("avgOffOccL322U", -1);
168  m_monObj->setVariable("avgOffOccL461U", -1);
169  m_monObj->setVariable("avgOffOccL462U", -1);
170  m_monObj->setVariable("avgOffOccL581U", -1);
171  m_monObj->setVariable("avgOffOccL582U", -1);
172  m_monObj->setVariable("avgOffOccL6101U", -1);
173  m_monObj->setVariable("avgOffOccL6102U", -1);
174  } else {
175  m_monObj->setVariable("avgOffOccL3U", avgOffOccL3UV[0]);
176  m_monObj->setVariable("avgOffOccL4U", avgOffOccL4UV[0]);
177  m_monObj->setVariable("avgOffOccL5U", avgOffOccL5UV[0]);
178  m_monObj->setVariable("avgOffOccL6U", avgOffOccL6UV[0]);
179  m_monObj->setVariable("avgOffOccL4X2U", avgOffOccL4X2UV[0]);
180  m_monObj->setVariable("avgOffOccL5X2U", avgOffOccL5X2UV[0]);
181  m_monObj->setVariable("avgOffOccL5X3U", avgOffOccL5X3UV[0]);
182  m_monObj->setVariable("avgOffOccL6X2U", avgOffOccL6X2UV[0]);
183  m_monObj->setVariable("avgOffOccL6X3U", avgOffOccL6X3UV[0]);
184  m_monObj->setVariable("avgOffOccL6X4U", avgOffOccL6X4UV[0]);
185  m_monObj->setVariable("avgOffOccL311U", avgOffOccL311UV[0]);
186  m_monObj->setVariable("avgOffOccL312U", avgOffOccL312UV[0]);
187  m_monObj->setVariable("avgOffOccL321U", avgOffOccL321UV[0]);
188  m_monObj->setVariable("avgOffOccL322U", avgOffOccL322UV[0]);
189  m_monObj->setVariable("avgOffOccL461U", avgOffOccL461UV[0]);
190  m_monObj->setVariable("avgOffOccL462U", avgOffOccL462UV[0]);
191  m_monObj->setVariable("avgOffOccL581U", avgOffOccL581UV[0]);
192  m_monObj->setVariable("avgOffOccL582U", avgOffOccL582UV[0]);
193  m_monObj->setVariable("avgOffOccL6101U", avgOffOccL6101UV[0]);
194  m_monObj->setVariable("avgOffOccL6102U", avgOffOccL6102UV[0]);
195  }
196 
197  if (h_zs5countsV == NULL || h_events == NULL) {
198  B2INFO("Histograms needed for Average Offline Occupancy on V side are not found");
199  m_monObj->setVariable("avgOffOccL3V", -1);
200  m_monObj->setVariable("avgOffOccL4V", -1);
201  m_monObj->setVariable("avgOffOccL5V", -1);
202  m_monObj->setVariable("avgOffOccL6V", -1);
203  m_monObj->setVariable("avgOffOccL4X2V", -1);
204  m_monObj->setVariable("avgOffOccL5X2V", -1);
205  m_monObj->setVariable("avgOffOccL5X3V", -1);
206  m_monObj->setVariable("avgOffOccL6X2V", -1);
207  m_monObj->setVariable("avgOffOccL6X3V", -1);
208  m_monObj->setVariable("avgOffOccL6X4V", -1);
209  m_monObj->setVariable("avgOffOccL311V", -1);
210  m_monObj->setVariable("avgOffOccL312V", -1);
211  m_monObj->setVariable("avgOffOccL321V", -1);
212  m_monObj->setVariable("avgOffOccL322V", -1);
213  m_monObj->setVariable("avgOffOccL461V", -1);
214  m_monObj->setVariable("avgOffOccL462V", -1);
215  m_monObj->setVariable("avgOffOccL581V", -1);
216  m_monObj->setVariable("avgOffOccL582V", -1);
217  m_monObj->setVariable("avgOffOccL6101V", -1);
218  m_monObj->setVariable("avgOffOccL6102V", -1);
219  } else {
220  m_monObj->setVariable("avgOffOccL3V", avgOffOccL3UV[1]);
221  m_monObj->setVariable("avgOffOccL4V", avgOffOccL4UV[1]);
222  m_monObj->setVariable("avgOffOccL5V", avgOffOccL5UV[1]);
223  m_monObj->setVariable("avgOffOccL6V", avgOffOccL6UV[1]);
224  m_monObj->setVariable("avgOffOccL4X2V", avgOffOccL4X2UV[1]);
225  m_monObj->setVariable("avgOffOccL5X2V", avgOffOccL5X2UV[1]);
226  m_monObj->setVariable("avgOffOccL5X3V", avgOffOccL5X3UV[1]);
227  m_monObj->setVariable("avgOffOccL6X2V", avgOffOccL6X2UV[1]);
228  m_monObj->setVariable("avgOffOccL6X3V", avgOffOccL6X3UV[1]);
229  m_monObj->setVariable("avgOffOccL6X4V", avgOffOccL6X4UV[1]);
230  m_monObj->setVariable("avgOffOccL311V", avgOffOccL311UV[1]);
231  m_monObj->setVariable("avgOffOccL312V", avgOffOccL312UV[1]);
232  m_monObj->setVariable("avgOffOccL321V", avgOffOccL321UV[1]);
233  m_monObj->setVariable("avgOffOccL322V", avgOffOccL322UV[1]);
234  m_monObj->setVariable("avgOffOccL461V", avgOffOccL461UV[1]);
235  m_monObj->setVariable("avgOffOccL462V", avgOffOccL462UV[1]);
236  m_monObj->setVariable("avgOffOccL581V", avgOffOccL581UV[1]);
237  m_monObj->setVariable("avgOffOccL582V", avgOffOccL582UV[1]);
238  m_monObj->setVariable("avgOffOccL6101V", avgOffOccL6101UV[1]);
239  m_monObj->setVariable("avgOffOccL6102V", avgOffOccL6102UV[1]);
240  }
241 
242 
243  // efficiency of cluster recontruction for U and V side
244  TH2F* h_found_tracksU = (TH2F*)findHist("SVDEfficiency/TrackHitsU");
245  TH2F* h_matched_clusU = (TH2F*)findHist("SVDEfficiency/MatchedHitsU");
246  TH2F* h_found_tracksV = (TH2F*)findHist("SVDEfficiency/TrackHitsV");
247  TH2F* h_matched_clusV = (TH2F*)findHist("SVDEfficiency/MatchedHitsV");
248 
249  m_c_avgEfficiency->Clear();
250  m_c_avgEfficiency->Divide(2, 2);
251  m_c_avgEfficiency->cd(1);
252  if (h_found_tracksU) h_found_tracksU->Draw("colz");
253  m_c_avgEfficiency->cd(2);
254  if (h_found_tracksV) h_found_tracksV->Draw("colz");
255  m_c_avgEfficiency->cd(3);
256  if (h_matched_clusU) h_matched_clusU->Draw("colz");
257  m_c_avgEfficiency->cd(4);
258  if (h_matched_clusV) h_matched_clusV->Draw("colz");
259 
260  // average efficiency in each layer for both side (U, V)
261  std::vector<float> avgEffL3 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 7, 2, 3);
262 
263  std::vector<float> avgEffL4 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 10, 5, 7);
264 
265  std::vector<float> avgEffL5 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 12, 9, 12);
266 
267  std::vector<float> avgEffL6 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 16, 14, 18);
268 
269  // average efficiency
270  std::vector<float> avgEffL3456 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 16, 2, 18);
271 
272  // average efficiency for mid plane +x L3.1.1, L3.1.2
273  std::vector<float> avgEffL31X = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 1, 2, 3);
274 
275  // average efficiency for mid plane +x L3.2.1, L3.2.2
276  std::vector<float> avgEffL32X = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 2, 2, 2, 3);
277 
278  // average efficiency for mid plane: L4.X.2
279  std::vector<float> avgEffL4X2 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 10, 6, 6);
280 
281  // average efficiency for mid plane: L5.X.2
282  std::vector<float> avgEffL5X2 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 12, 10, 10);
283 
284  // average efficiency for mid plane: L5.X.3
285  std::vector<float> avgEffL5X3 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 12, 11, 11);
286 
287  // average efficiency for mid plane: L6.X.2
288  std::vector<float> avgEffL6X2 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 16, 15, 15);
289 
290  // average efficiency for mid plane: L6.X.3
291  std::vector<float> avgEffL6X3 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 16, 16, 16);
292 
293  // average efficiency for mid plane: L6.X.4
294  std::vector<float> avgEffL6X4 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 16, 17, 17);
295 
296  // average efficiency for high occupancy sensors
297 
298  std::vector<float> avgEffL311UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 1, 2,
299  2); // L3.1.1
300 
301  std::vector<float> avgEffL312UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 1, 1, 3,
302  3); // L3.1.2
303 
304  std::vector<float> avgEffL321UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 2, 2, 2,
305  2); // L3.2.1
306 
307  std::vector<float> avgEffL322UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 2, 2, 3,
308  3); // L3.2.2
309 
310  std::vector<float> avgEffL461UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 6, 6, 5,
311  5); // L4.6.1
312 
313  std::vector<float> avgEffL462UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 6, 6, 6,
314  6); // L4.6.2
315 
316  std::vector<float> avgEffL581UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 8, 8, 9,
317  9); // L5.8.1
318 
319  std::vector<float> avgEffL582UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 8, 8, 10,
320  10); // L5.8.2
321 
322  std::vector<float> avgEffL6101UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 10, 10, 14,
323  14); // L6.10.1
324 
325  std::vector<float> avgEffL6102UV = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 10, 10, 15,
326  15); // L6.10.2
327 
328  // setting monitoring variables
329  if (h_matched_clusU == NULL || h_found_tracksU == NULL) {
330  B2INFO("Histograms needed for Average Efficiency on U side are not found");
331  m_monObj->setVariable("avgEffL3U", -1);
332  m_monObj->setVariable("avgEffL4U", -1);
333  m_monObj->setVariable("avgEffL5U", -1);
334  m_monObj->setVariable("avgEffL6U", -1);
335  m_monObj->setVariable("avgEffL3456U", -1);
336  m_monObj->setVariable("avgEffL31XU", -1);
337  m_monObj->setVariable("avgEffL32XU", -1);
338  m_monObj->setVariable("avgEffL4X2U", -1);
339  m_monObj->setVariable("avgEffL5X2U", -1);
340  m_monObj->setVariable("avgEffL5X3U", -1);
341  m_monObj->setVariable("avgEffL6X2U", -1);
342  m_monObj->setVariable("avgEffL6X3U", -1);
343  m_monObj->setVariable("avgEffL6X4U", -1);
344  m_monObj->setVariable("avgEffL311U", -1);
345  m_monObj->setVariable("avgEffL312U", -1);
346  m_monObj->setVariable("avgEffL321U", -1);
347  m_monObj->setVariable("avgEffL322U", -1);
348  m_monObj->setVariable("avgEffL461U", -1);
349  m_monObj->setVariable("avgEffL462U", -1);
350  m_monObj->setVariable("avgEffL581U", -1);
351  m_monObj->setVariable("avgEffL582U", -1);
352  m_monObj->setVariable("avgEffL6101U", -1);
353  m_monObj->setVariable("avgEffL6102U", -1);
354  } else {
355  m_monObj->setVariable("avgEffL3U", avgEffL3[0]);
356  m_monObj->setVariable("avgEffL4U", avgEffL4[0]);
357  m_monObj->setVariable("avgEffL5U", avgEffL5[0]);
358  m_monObj->setVariable("avgEffL6U", avgEffL6[0]);
359  m_monObj->setVariable("avgEffL3456U", avgEffL3456[0]);
360  m_monObj->setVariable("avgEffL31XU", avgEffL31X[0]);
361  m_monObj->setVariable("avgEffL32XU", avgEffL32X[0]);
362  m_monObj->setVariable("avgEffL4X2U", avgEffL4X2[0]);
363  m_monObj->setVariable("avgEffL5X2U", avgEffL5X2[0]);
364  m_monObj->setVariable("avgEffL5X3U", avgEffL5X3[0]);
365  m_monObj->setVariable("avgEffL6X2U", avgEffL6X2[0]);
366  m_monObj->setVariable("avgEffL6X3U", avgEffL6X3[0]);
367  m_monObj->setVariable("avgEffL6X4U", avgEffL6X4[0]);
368  m_monObj->setVariable("avgEffL311U", avgEffL311UV[0]);
369  m_monObj->setVariable("avgEffL312U", avgEffL312UV[0]);
370  m_monObj->setVariable("avgEffL321U", avgEffL321UV[0]);
371  m_monObj->setVariable("avgEffL322U", avgEffL322UV[0]);
372  m_monObj->setVariable("avgEffL461U", avgEffL461UV[0]);
373  m_monObj->setVariable("avgEffL462U", avgEffL462UV[0]);
374  m_monObj->setVariable("avgEffL581U", avgEffL581UV[0]);
375  m_monObj->setVariable("avgEffL582U", avgEffL582UV[0]);
376  m_monObj->setVariable("avgEffL6101U", avgEffL6101UV[0]);
377  m_monObj->setVariable("avgEffL6102U", avgEffL6102UV[0]);
378  }
379 
380  if (h_matched_clusV == NULL || h_found_tracksV == NULL) {
381  B2INFO("Histograms needed for Average Efficiency on V side are not found");
382  m_monObj->setVariable("avgEffL3V", -1);
383  m_monObj->setVariable("avgEffL4V", -1);
384  m_monObj->setVariable("avgEffL5V", -1);
385  m_monObj->setVariable("avgEffL6V", -1);
386  m_monObj->setVariable("avgEffL3456V", -1);
387  m_monObj->setVariable("avgEffL31XV", -1);
388  m_monObj->setVariable("avgEffL32XV", -1);
389  m_monObj->setVariable("avgEffL4X2V", -1);
390  m_monObj->setVariable("avgEffL5X2V", -1);
391  m_monObj->setVariable("avgEffL5X3V", -1);
392  m_monObj->setVariable("avgEffL6X2V", -1);
393  m_monObj->setVariable("avgEffL6X3V", -1);
394  m_monObj->setVariable("avgEffL6X4V", -1);
395  m_monObj->setVariable("avgEffL311V", -1);
396  m_monObj->setVariable("avgEffL312V", -1);
397  m_monObj->setVariable("avgEffL321V", -1);
398  m_monObj->setVariable("avgEffL322V", -1);
399  m_monObj->setVariable("avgEffL461V", -1);
400  m_monObj->setVariable("avgEffL462V", -1);
401  m_monObj->setVariable("avgEffL581V", -1);
402  m_monObj->setVariable("avgEffL582V", -1);
403  m_monObj->setVariable("avgEffL6101V", -1);
404  m_monObj->setVariable("avgEffL6102V", -1);
405  } else {
406  m_monObj->setVariable("avgEffL3V", avgEffL3[1]);
407  m_monObj->setVariable("avgEffL4V", avgEffL4[1]);
408  m_monObj->setVariable("avgEffL5V", avgEffL5[1]);
409  m_monObj->setVariable("avgEffL6V", avgEffL6[1]);
410  m_monObj->setVariable("avgEffL3456V", avgEffL3456[1]);
411  m_monObj->setVariable("avgEffL31XV", avgEffL31X[1]);
412  m_monObj->setVariable("avgEffL32XV", avgEffL32X[1]);
413  m_monObj->setVariable("avgEffL4X2V", avgEffL4X2[1]);
414  m_monObj->setVariable("avgEffL5X2V", avgEffL5X2[1]);
415  m_monObj->setVariable("avgEffL5X3V", avgEffL5X3[1]);
416  m_monObj->setVariable("avgEffL6X2V", avgEffL6X2[1]);
417  m_monObj->setVariable("avgEffL6X3V", avgEffL6X3[1]);
418  m_monObj->setVariable("avgEffL6X4V", avgEffL6X4[1]);
419  m_monObj->setVariable("avgEffL311U", avgEffL311UV[1]);
420  m_monObj->setVariable("avgEffL312U", avgEffL312UV[1]);
421  m_monObj->setVariable("avgEffL321U", avgEffL321UV[1]);
422  m_monObj->setVariable("avgEffL322U", avgEffL322UV[1]);
423  m_monObj->setVariable("avgEffL461U", avgEffL461UV[1]);
424  m_monObj->setVariable("avgEffL462U", avgEffL462UV[1]);
425  m_monObj->setVariable("avgEffL581U", avgEffL581UV[1]);
426  m_monObj->setVariable("avgEffL582U", avgEffL582UV[1]);
427  m_monObj->setVariable("avgEffL6101U", avgEffL6101UV[1]);
428  m_monObj->setVariable("avgEffL6102U", avgEffL6102UV[1]);
429  }
430 
431 
432  // MPV cluster charge for clusters on track
433  TH1F* h_clusterCharge_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeU3");
434  TH1F* h_clusterCharge_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeV3");
435  TH1F* h_clusterCharge_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeU456");
436  TH1F* h_clusterCharge_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeV456");
437 
438  m_c_MPVChargeClusterOnTrack->Clear();
439  m_c_MPVChargeClusterOnTrack->Divide(2, 2);
440  m_c_MPVChargeClusterOnTrack->cd(1);
441  if (h_clusterCharge_L3U) h_clusterCharge_L3U->Draw();
442  m_c_MPVChargeClusterOnTrack->cd(2);
443  if (h_clusterCharge_L3V) h_clusterCharge_L3V->Draw();
444  m_c_MPVChargeClusterOnTrack->cd(3);
445  if (h_clusterCharge_L456U) h_clusterCharge_L456U->Draw();
446  m_c_MPVChargeClusterOnTrack->cd(4);
447  if (h_clusterCharge_L456U) h_clusterCharge_L456U->Draw();
448 
449  // find abscissa of max Y in histograms
450  float MPVClusterChargeL3U = xForMaxY(h_clusterCharge_L3U);
451  float MPVClusterChargeL3V = xForMaxY(h_clusterCharge_L3V);
452  float MPVClusterChargeL456U = xForMaxY(h_clusterCharge_L456U);
453  float MPVClusterChargeL456V = xForMaxY(h_clusterCharge_L456V);
454 
455  if (h_clusterCharge_L3U == NULL || h_clusterCharge_L456U == NULL) {
456  B2INFO("Histograms needed for MPV cluster charge on U side are not found");
457  m_monObj->setVariable("MPVClusterChargeL3U", -1);
458  m_monObj->setVariable("MPVClusterChargeL456U", -1);
459  } else {
460  m_monObj->setVariable("MPVClusterChargeL3U", MPVClusterChargeL3U);
461  m_monObj->setVariable("MPVClusterChargeL456U", MPVClusterChargeL456U);
462  }
463 
464  if (h_clusterCharge_L3V == NULL || h_clusterCharge_L456V == NULL) {
465  B2INFO("Histograms needed for MPV cluster charge on V side are not found");
466  m_monObj->setVariable("MPVClusterChargeL3V", -1);
467  m_monObj->setVariable("MPVClusterChargeL456V", -1);
468  } else {
469  m_monObj->setVariable("MPVClusterChargeL3V", MPVClusterChargeL3V);
470  m_monObj->setVariable("MPVClusterChargeL456V", MPVClusterChargeL456V);
471  }
472 
473 
474  // MPV SNR for the clusters on track
475  TH1F* h_clusterSNR_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRU3");
476  TH1F* h_clusterSNR_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRV3");
477  TH1F* h_clusterSNR_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRU456");
478  TH1F* h_clusterSNR_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRV456");
479 
480  m_c_MPVSNRClusterOnTrack->Clear();
481  m_c_MPVSNRClusterOnTrack->Divide(2, 2);
482  m_c_MPVSNRClusterOnTrack->cd(1);
483  if (h_clusterSNR_L3U) h_clusterSNR_L3U->Draw();
484  m_c_MPVSNRClusterOnTrack->cd(2);
485  if (h_clusterSNR_L3V) h_clusterSNR_L3V->Draw();
486  m_c_MPVSNRClusterOnTrack->cd(3);
487  if (h_clusterSNR_L456U) h_clusterSNR_L456U->Draw();
488  m_c_MPVSNRClusterOnTrack->cd(4);
489  if (h_clusterSNR_L456V) h_clusterSNR_L456V->Draw();
490 
491  float MPVClusterSNRL3U = xForMaxY(h_clusterSNR_L3U);
492  float MPVClusterSNRL3V = xForMaxY(h_clusterSNR_L3V);
493  float MPVClusterSNRL456U = xForMaxY(h_clusterSNR_L456U);
494  float MPVClusterSNRL456V = xForMaxY(h_clusterSNR_L456V);
495 
496  if (h_clusterSNR_L3U == NULL || h_clusterSNR_L456U == NULL) {
497  B2INFO("Histograms needed for MPV cluster SNR on U side are not found");
498  m_monObj->setVariable("MPVClusterSNRL3U", -1);
499  m_monObj->setVariable("MPVClusterSNRL456U", -1);
500  } else {
501  m_monObj->setVariable("MPVClusterSNRL3U", MPVClusterSNRL3U);
502  m_monObj->setVariable("MPVClusterSNRL456U", MPVClusterSNRL456U);
503  }
504 
505  if (h_clusterSNR_L3V == NULL || h_clusterSNR_L456V == NULL) {
506  B2INFO("Histograms needed for MPV cluster SNR on V side are not found");
507  m_monObj->setVariable("MPVClusterSNRL3V", -1);
508  m_monObj->setVariable("MPVClusterSNRL456V", -1);
509  } else {
510  m_monObj->setVariable("MPVClusterSNRL3V", MPVClusterSNRL3V);
511  m_monObj->setVariable("MPVClusterSNRL456V", MPVClusterSNRL456V);
512  }
513 
514 
515  // MPV SVD cluster time for the clusters on track
516  TH1F* h_clusterTime_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeU3");
517  TH1F* h_clusterTime_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeV3");
518  TH1F* h_clusterTime_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeU456");
519  TH1F* h_clusterTime_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeV456");
520 
521  m_c_MPVTimeClusterOnTrack->Clear();
522  m_c_MPVTimeClusterOnTrack->Divide(2, 2);
523  m_c_MPVTimeClusterOnTrack->cd(1);
524  if (h_clusterTime_L3U) h_clusterTime_L3U->Draw();
525  m_c_MPVTimeClusterOnTrack->cd(2);
526  if (h_clusterTime_L3V) h_clusterTime_L3V->Draw();
527  m_c_MPVTimeClusterOnTrack->cd(3);
528  if (h_clusterTime_L456U) h_clusterTime_L456U->Draw();
529  m_c_MPVTimeClusterOnTrack->cd(4);
530  if (h_clusterTime_L456V) h_clusterTime_L456V->Draw();
531 
532  float MPVClusterTimeL3U = xForMaxY(h_clusterTime_L3U);
533  float MPVClusterTimeL3V = xForMaxY(h_clusterTime_L3V);
534  float MPVClusterTimeL456U = xForMaxY(h_clusterTime_L456U);
535  float MPVClusterTimeL456V = xForMaxY(h_clusterTime_L456V);
536  float FWHMClusterTimeL3U = histFWHM(h_clusterTime_L3U);
537  float FWHMClusterTimeL3V = histFWHM(h_clusterTime_L3V);
538  float FWHMClusterTimeL456U = histFWHM(h_clusterTime_L456U);
539  float FWHMClusterTimeL456V = histFWHM(h_clusterTime_L456V);
540 
541  if (h_clusterTime_L3U == NULL || h_clusterTime_L456U == NULL) {
542  B2INFO("Histograms needed for MPV cluster time on U side are not found");
543  m_monObj->setVariable("MPVClusterTimeL3U", -1);
544  m_monObj->setVariable("MPVClusterTimeL456U", -1);
545  m_monObj->setVariable("FWHMClusterTimeL3U", -1);
546  m_monObj->setVariable("FWHMClusterTimeL456U", -1);
547  } else {
548  m_monObj->setVariable("MPVClusterTimeL3U", MPVClusterTimeL3U);
549  m_monObj->setVariable("MPVClusterTimeL456U", MPVClusterTimeL456U);
550  m_monObj->setVariable("FWHMClusterTimeL3U", FWHMClusterTimeL3U);
551  m_monObj->setVariable("FWHMClusterTimeL456U", FWHMClusterTimeL456U);
552  }
553 
554  if (h_clusterTime_L3V == NULL || h_clusterTime_L456V == NULL) {
555  B2INFO("Histograms needed for MPV cluster time on V side are not found");
556  m_monObj->setVariable("MPVClusterTimeL3V", -1);
557  m_monObj->setVariable("MPVClusterTimeL456V", -1);
558  m_monObj->setVariable("FWHMClusterTimeL3V", -1);
559  m_monObj->setVariable("FWHMClusterTimeL456V", -1);
560  } else {
561  m_monObj->setVariable("MPVClusterTimeL3V", MPVClusterTimeL3V);
562  m_monObj->setVariable("MPVClusterTimeL456V", MPVClusterTimeL456V);
563  m_monObj->setVariable("FWHMClusterTimeL3V", FWHMClusterTimeL3V);
564  m_monObj->setVariable("FWHMClusterTimeL456V", FWHMClusterTimeL456V);
565  }
566 
567 
568  // average maxBin for clusters on track
569  TH1F* h_maxBinU = (TH1F*)findHist("SVDClsTrk/SVDTRK_StripMaxBinUAll");
570  TH1F* h_maxBinV = (TH1F*)findHist("SVDClsTrk/SVDTRK_StripMaxBinVAll");
571 
572  m_c_avgMaxBinClusterOnTrack->Clear();
573  m_c_avgMaxBinClusterOnTrack->Divide(2, 1);
574  m_c_avgMaxBinClusterOnTrack->cd(1);
575  if (h_maxBinU) h_maxBinU->Draw();
576  m_c_avgMaxBinClusterOnTrack->cd(2);
577  if (h_maxBinV) h_maxBinV->Draw();
578 
579  if (h_maxBinU == NULL) {
580  B2INFO("Histogram needed for Average MaxBin on U side is not found");
581  m_monObj->setVariable("avgMaxBinU", -1);
582  } else {
583  float avgMaxBinU = h_maxBinU->GetMean();
584  m_monObj->setVariable("avgMaxBinU", avgMaxBinU);
585  }
586 
587  if (h_maxBinV == NULL) {
588  B2INFO("Histogram needed for Average MaxBin on V side is not found");
589  m_monObj->setVariable("avgMaxBinV", -1);
590  } else {
591  float avgMaxBinV = h_maxBinV->GetMean();
592  m_monObj->setVariable("avgMaxBinV", avgMaxBinV);
593  }
594 
595  B2INFO("DQMHistAnalysisSVDGeneral: endRun called");
596 }
597 
598 
599 void DQMHistAnalysisSVDOnMiraBelleModule::terminate()
600 {
601  B2INFO("DQMHistAnalysisSVDOnMiraBelle: terminate called");
602 }
603 
604 
605 std::vector<float> DQMHistAnalysisSVDOnMiraBelleModule::highOccupancySensor(int iLayer, TH1F* hU, TH1F* hV, int iBin,
606  int nEvents) const
607 {
608  int nStripsV = -1;
609  if (iLayer == 3) {
610  nStripsV = 768;
611  } else if (iLayer >= 4 && iLayer <= 6) {
612  nStripsV = 512;
613  } else {
614  B2DEBUG(20, "Layer out of range [3,6].");
615  }
616  std::vector<float> avgOffOccUV(2, 0.0);
617  avgOffOccUV[0] = hU->GetBinContent(iBin) / 768 / nEvents * 100;
618  avgOffOccUV[1] = hV->GetBinContent(iBin) / nStripsV / nEvents * 100;
619  return avgOffOccUV;
620 }
621 
622 
623 std::vector<float> DQMHistAnalysisSVDOnMiraBelleModule::avgOccupancyUV(int iLayer, TH1F* hU, TH1F* hV, int min, int max, int offset,
624  int step, int nEvents) const
625 {
626  int nStripsV = -1;
627  if (iLayer == 3) {
628  nStripsV = 768;
629  } else if (iLayer >= 4 && iLayer <= 6) {
630  nStripsV = 512;
631  } else {
632  B2DEBUG(20, "Layer out of range [3,6].");
633  }
634  std::vector<float> avgOffOccUV(2, 0.0);
635  for (int bin = min; bin < max; bin++) {
636  avgOffOccUV[0] += hU->GetBinContent(offset + step * bin) / 768 * 100;
637  avgOffOccUV[1] += hV->GetBinContent(offset + step * bin) / nStripsV * 100;
638  }
639  avgOffOccUV[0] /= ((max - min) * nEvents);
640  avgOffOccUV[1] /= ((max - min) * nEvents);
641  return avgOffOccUV;
642 }
643 
644 
645 std::vector<float> DQMHistAnalysisSVDOnMiraBelleModule::avgEfficiencyUV(TH2F* hMCU, TH2F* hMCV, TH2F* hFTU, TH2F* hFTV, int minX,
646  int maxX, int minY, int maxY) const
647 {
648  std::vector<float> avgEffUV(2, 0.0);
649  std::vector<float> sumMatchedClustersUV(2, 0.0);
650  std::vector<float> sumFoundTracksUV(2, 0.0);
651  for (int binX = minX; binX < maxX + 1; binX++) {
652  for (int binY = minY; binY < maxY + 1; binY++) {
653  int binXY = hMCV->GetBin(binX, binY);
654  sumMatchedClustersUV[0] += hMCU->GetBinContent(binXY);
655  sumMatchedClustersUV[1] += hMCV->GetBinContent(binXY);
656  sumFoundTracksUV[0] += hFTU->GetBinContent(binXY);
657  sumFoundTracksUV[1] += hFTV->GetBinContent(binXY);
658  }
659  }
660  if (sumFoundTracksUV[0] > 0) {
661  avgEffUV[0] = sumMatchedClustersUV[0] / sumFoundTracksUV[0] * 100;
662  } else {
663  avgEffUV[0] = -1;
664  }
665  if (sumFoundTracksUV[1] > 0) {
666  avgEffUV[1] = sumMatchedClustersUV[1] / sumFoundTracksUV[1] * 100;
667  } else {
668  avgEffUV[1] = -1;
669  }
670  return avgEffUV;
671 }
672 
673 float DQMHistAnalysisSVDOnMiraBelleModule::xForMaxY(TH1F* h) const
674 {
675  int maxY = h->GetMaximumBin();
676  float xMaxY = h->GetXaxis()->GetBinCenter(maxY);
677  return xMaxY;
678 }
679 
680 float DQMHistAnalysisSVDOnMiraBelleModule::histFWHM(TH1F* h) const
681 {
682  int bin1 = h->FindFirstBinAbove(h->GetMaximum() / 2);
683  int bin2 = h->FindLastBinAbove(h->GetMaximum() / 2);
684  float fwhm = h->GetBinCenter(bin2) - h->GetBinCenter(bin1);
685  return fwhm;
686 }
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisSVDOnMiraBelleModule
Class derived from HistoModule, for SVD monitoring variables at MiraBelle.
Definition: DQMHistAnalysisSVDOnMiraBelle.h:37
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27