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