Belle II Software development
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 <numeric>
10#include <limits>
11#include <vxd/geometry/GeoCache.h>
12#include <vxd/dataobjects/VxdID.h>
13#include <svd/geometry/SensorInfo.h>
14#include <vxd/geometry/SensorInfoBase.h>
15#include <vxd/geometry/GeoTools.h>
16#include <framework/datastore/StoreObjPtr.h>
17#include <framework/datastore/StoreArray.h>
18
19
20#include <dqm/analysis/modules/DQMHistAnalysisSVDOnMiraBelle.h>
21
22using namespace std;
23using namespace Belle2;
24
25//-----------------------------------------------------------------
26// Register the Module
27//-----------------------------------------------------------------
28REG_MODULE(DQMHistAnalysisSVDOnMiraBelle);
29
30//-----------------------------------------------------------------
31// Implementation
32//-----------------------------------------------------------------
33
36{
37 setDescription("DQM Analysis Module that extracts monitoring variables from SVD DQM histograms and provides input to MiraBelle.");
39 B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: Constructor done.");
40}
41
43
45{
46 gROOT->cd();
47
49
50 // add MonitoringObject
52
53 // list of canvases to be added to MonitoringObject
54 m_c_avgEfficiency = new TCanvas("svd_avgEfficiency", "matched clusters and found tracks", 0, 0, 800, 600);
55 m_c_avgOffOccupancy = new TCanvas("svd_avgOffOccupancy", "strips", 0, 0, 800, 600);
56 m_c_MPVChargeClusterOnTrack = new TCanvas("svd_MPVChargeClusterOnTrack", "charge from Clusters on Track Charge", 0, 0, 400, 400);
57 m_c_MPVSNRClusterOnTrack = new TCanvas("svd_MPVSNRClusterOnTrack", "SNR from Clusters on Track Charge", 0, 0, 400, 400);
58 m_c_MPVTimeClusterOnTrack = new TCanvas("svd_MPVTimeClusterOnTrack", "time from Clusters on Track Charge", 0, 0, 400, 400);
59 m_c_avgMaxBinClusterOnTrack = new TCanvas("svd_avgMaxBin", "average MaxBin", 0, 0, 800, 600);
60 m_c_MeanSVDEventT0 = new TCanvas("svd_MeanSVDEventT0", "Mean Event T0 from SVD for all samples", 0, 0, 400, 400);
61
62 // add canvases used to create monitoring variables to MonitoringObject
63 m_monObj->addCanvas(m_c_avgEfficiency);
64 m_monObj->addCanvas(m_c_avgOffOccupancy);
69 m_monObj->addCanvas(m_c_MeanSVDEventT0);
70
72
73 //collect the list of all SVD Modules in the geometry here
74 std::vector<VxdID> sensors = geo.getListOfSensors();
75 for (VxdID& aVxdID : sensors) {
76 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
77 // B2INFO("VXD " << aVxdID);
78 if (info.getType() != VXD::SensorInfoBase::SVD) continue;
79 m_SVDModules.push_back(aVxdID); // reorder, sort would be better
80 }
81 std::sort(m_SVDModules.begin(), m_SVDModules.end()); // back to natural order
82
83 if (!m_svdPlotsConfig.isValid())
84 B2FATAL("no valid configuration found for SVD reconstruction");
85 else {
86 B2DEBUG(20, "SVDRecoConfiguration: from now on we are using " << m_svdPlotsConfig->get_uniqueID());
87 //read back from payload
88 m_listOfSensorsToMonitor = m_svdPlotsConfig->getListOfSensors();
89 }
90
91 B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: initialized.");
92}
93
95{
96 B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: beginRun called.");
97}
98
100{
101 B2DEBUG(20, "DQMHistAnalysisSVDOnMiraBelle: event called.");
102}
103
105{
106 float nan = numeric_limits<float>::quiet_NaN();
107
108 // ladder label
109 std::vector<string> ladderLabel = {"L3.X.1", "L3.X.2", "L4.X.1", "L4.X.2", "L4.X.3", "L5.X.1", "L5.X.2", "L5.X.3", "L5.X.4",
110 "L6.X.1", "L6.X.2", "L6.X.3", "L6.X.4", "L6.X.5"
111 };
112
113 // sensors to monitored
114 // for high occupancy sensors
115 // and for low DCDC sensors
116 std::vector<string> sensorLabel = {"L3.1.1", "L3.1.2", "L3.2.1", "L3.2.2", "L4.6.1", "L4.6.2", "L5.8.1", "L5.8.2", "L6.10.1", "L6.10.2",
117 "L4.1.1", "L4.10.2", "L4.3.3", "L5.1.3", "L5.1.4", "L5.9.2", "L5.9.4", "L6.4.3", "L6.6.4", "L6.10.3", "L6.11.5", "L6.12.4"
118 };
119
120 // offline occupancy - integrated number of ZS5 fired strips
121 TH1F* h_zs5countsU = (TH1F*)findHist("SVDExpReco/SVDDQM_StripCountsU"); // made by SVDDQMExperssRecoModule
122 TH1F* h_zs5countsV = (TH1F*)findHist("SVDExpReco/SVDDQM_StripCountsV");
123 TH1F* h_events = (TH1F*)findHist("SVDExpReco/SVDDQM_nEvents");
124
125 // adding histograms to canvas
126 m_c_avgOffOccupancy->Clear();
127 m_c_avgOffOccupancy->Divide(2, 2);
128 m_c_avgOffOccupancy->cd(1);
129 if (h_zs5countsU) h_zs5countsU->Draw("colz");
130 m_c_avgOffOccupancy->cd(2);
131 if (h_zs5countsV) h_zs5countsV->Draw("colz");
132 m_c_avgOffOccupancy->cd(3);
133 if (h_events) h_events->Draw("colz");
134
135 int nE = 0;
136 if (h_events) nE = h_events->GetEntries(); // number of events for all "clusters on track" histograms
137
138 // setting monitoring variables
139 if (h_zs5countsU == NULL || h_zs5countsV == NULL || h_events == NULL) {
140 if (h_zs5countsU == NULL) {
141 B2INFO("Histograms needed for Average Offline Occupancy on U side are not found");
142 }
143 if (h_zs5countsV == NULL) {
144 B2INFO("Histograms needed for Average Offline Occupancy on V side are not found");
145 }
146 } else {
147 // average occupancy for each layer
148 std::pair<float, float> avgOffOccL3UV = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, 3);
149 SetVariable(avgOffOccL3UV);
150
151 std::pair<float, float> avgOffOccL4UV = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, 4);
152 SetVariable(avgOffOccL4UV);
153
154 std::pair<float, float> avgOffOccL5UV = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, 5);
155 SetVariable(avgOffOccL5UV);
156
157 std::pair<float, float> avgOffOccL6UV = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, 6);
158 SetVariable(avgOffOccL6UV);
159
160 // average occupancy for each layer for group Id0
161 std::pair<float, float> avgOffGrpId0OccL3UV = avgOccupancyGrpId0UV(3, nE);
162 SetVariable(avgOffGrpId0OccL3UV);
163
164 std::pair<float, float> avgOffGrpId0OccL4UV = avgOccupancyGrpId0UV(4, nE);
165 SetVariable(avgOffGrpId0OccL4UV);
166
167 std::pair<float, float> avgOffGrpId0OccL5UV = avgOccupancyGrpId0UV(5, nE);
168 SetVariable(avgOffGrpId0OccL5UV);
169
170 std::pair<float, float> avgOffGrpId0OccL6UV = avgOccupancyGrpId0UV(6, nE);
171 SetVariable(avgOffGrpId0OccL6UV);
172
173 // occupancy averaged over ladders
174 for (const auto& it : ladderLabel) {
175 string sensorDescr = it;
176 int layer = 0;
177 int sensor = 0;
178 sscanf(it.c_str(), "L%d.X.%d", &layer, &sensor);
179 std::pair<float, float> avgOffOccL = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, layer, -1, sensor);
180 addVariable(Form("avgOffOccL%dX%dUV", layer, sensor), avgOffOccL);
181 }
182
183 // average occupancy for high occupancy sensors
184 for (const auto& it : sensorLabel) {
185 string sensorDescr = it;
186 int layer = 0;
187 int ladder = 0;
188 int sensor = 0;
189 sscanf(it.c_str(), "L%d.%d.%d", &layer, &ladder, &sensor);
190 std::pair<float, float> avgOffOccL = avgOccupancyUV(h_zs5countsU, h_zs5countsV, nE, layer, ladder, sensor);
191 addVariable(Form("avgOffOccL%d%d%dUV", layer, ladder, sensor), avgOffOccL);
192 }
193 }
194
195
196 // efficiency of cluster reconstruction for U and V side
197 TH2F* h_found_tracksU = (TH2F*)findHist("SVDEfficiency/TrackHitsU");
198 TH2F* h_matched_clusU = (TH2F*)findHist("SVDEfficiency/MatchedHitsU");
199 TH2F* h_found_tracksV = (TH2F*)findHist("SVDEfficiency/TrackHitsV");
200 TH2F* h_matched_clusV = (TH2F*)findHist("SVDEfficiency/MatchedHitsV");
201
202 m_c_avgEfficiency->Clear();
203 m_c_avgEfficiency->Divide(2, 2);
204 m_c_avgEfficiency->cd(1);
205 if (h_found_tracksU) h_found_tracksU->Draw("colz");
206 m_c_avgEfficiency->cd(2);
207 if (h_found_tracksV) h_found_tracksV->Draw("colz");
208 m_c_avgEfficiency->cd(3);
209 if (h_matched_clusU) h_matched_clusU->Draw("colz");
210 m_c_avgEfficiency->cd(4);
211 if (h_matched_clusV) h_matched_clusV->Draw("colz");
212
213 // setting monitoring variables
214 if (h_matched_clusU == NULL || h_matched_clusV == NULL || h_found_tracksU == NULL) {
215 if (h_matched_clusU == NULL) {
216 B2INFO("Histograms needed for Average Efficiency on U side are not found");
217 }
218 if (h_matched_clusV == NULL) {
219 B2INFO("Histograms needed for Average Efficiency on V side are not found");
220 }
221 } else {
222 // average efficiency in each layer for both side (U, V)
223 std::pair<float, float> avgEffL3 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 3);
224 SetVariable(avgEffL3);
225
226 std::pair<float, float> avgEffL4 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 4);
227 SetVariable(avgEffL4);
228
229 std::pair<float, float> avgEffL5 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 5);
230 SetVariable(avgEffL5);
231
232 std::pair<float, float> avgEffL6 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, 6);
233 SetVariable(avgEffL6);
234
235 // average efficiency for all layers
236 std::pair<float, float> avgEffL3456 = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV);
237 SetVariable(avgEffL3456);
238
239 // efficiency averaged over ladders
240 for (const auto& it : ladderLabel) {
241 string sensorDescr = it;
242 int layer = 0;
243 int sensor = 0;
244 sscanf(it.c_str(), "L%d.X.%d", &layer, &sensor);
245 std::pair<float, float> avgEffL = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, layer, -1,
246 sensor);
247 addVariable(Form("avgEffL%dX%dUV", layer, sensor), avgEffL);
248 }
249
250 // average efficiency for high occupancy sensors and
251 // average efficiency for low DCDC
252 for (const auto& it : sensorLabel) {
253 string sensorDescr = it;
254 int layer = 0;
255 int ladder = 0;
256 int sensor = 0;
257 sscanf(it.c_str(), "L%d.%d.%d", &layer, &ladder, &sensor);
258 std::pair<float, float> avgEffL = avgEfficiencyUV(h_matched_clusU, h_matched_clusV, h_found_tracksU, h_found_tracksV, layer, ladder,
259 sensor);
260 addVariable(Form("avgEffL%d%d%dUV", layer, ladder, sensor), avgEffL);
261 }
262 }
263
264 // MPV cluster charge for clusters on track
265 TH1F* h_clusterCharge_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeU3");
266 TH1F* h_clusterCharge_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeV3");
267 TH1F* h_clusterCharge_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeU456");
268 TH1F* h_clusterCharge_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterChargeV456");
269
271 m_c_MPVChargeClusterOnTrack->Divide(2, 2);
273 if (h_clusterCharge_L3U) h_clusterCharge_L3U->Draw();
275 if (h_clusterCharge_L3V) h_clusterCharge_L3V->Draw();
277 if (h_clusterCharge_L456U) h_clusterCharge_L456U->Draw();
279 if (h_clusterCharge_L456U) h_clusterCharge_L456U->Draw();
280
281 // find abscissa of max Y in histograms
282 float MPVClusterChargeL3U = nan;
283 if (h_clusterCharge_L3U)
284 if (h_clusterCharge_L3U->GetEntries() != 0)
285 MPVClusterChargeL3U = xForMaxY(h_clusterCharge_L3U);
286 float MPVClusterChargeL3V = nan;
287 if (h_clusterCharge_L3V)
288 if (h_clusterCharge_L3V->GetEntries() != 0)
289 MPVClusterChargeL3V = xForMaxY(h_clusterCharge_L3V);
290 float MPVClusterChargeL456U = nan;
291 if (h_clusterCharge_L456U)
292 if (h_clusterCharge_L456U->GetEntries() != 0)
293 MPVClusterChargeL456U = xForMaxY(h_clusterCharge_L456U);
294 float MPVClusterChargeL456V = nan;
295 if (h_clusterCharge_L456V)
296 if (h_clusterCharge_L456V->GetEntries() != 0)
297 MPVClusterChargeL456V = xForMaxY(h_clusterCharge_L456V);
298
299 if (h_clusterCharge_L3U == NULL || h_clusterCharge_L456U == NULL) {
300 B2INFO("Histograms needed for MPV cluster charge on U side are not found");
301 } else {
302 m_monObj->setVariable("MPVClusterChargeL3U", MPVClusterChargeL3U);
303 m_monObj->setVariable("MPVClusterChargeL456U", MPVClusterChargeL456U);
304 }
305
306 if (h_clusterCharge_L3V == NULL || h_clusterCharge_L456V == NULL) {
307 B2INFO("Histograms needed for MPV cluster charge on V side are not found");
308 } else {
309 m_monObj->setVariable("MPVClusterChargeL3V", MPVClusterChargeL3V);
310 m_monObj->setVariable("MPVClusterChargeL456V", MPVClusterChargeL456V);
311 }
312
313
314 // MPV SNR for the clusters on track
315 TH1F* h_clusterSNR_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRU3");
316 TH1F* h_clusterSNR_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRV3");
317 TH1F* h_clusterSNR_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRU456");
318 TH1F* h_clusterSNR_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterSNRV456");
319
321 m_c_MPVSNRClusterOnTrack->Divide(2, 2);
323 if (h_clusterSNR_L3U) h_clusterSNR_L3U->Draw();
325 if (h_clusterSNR_L3V) h_clusterSNR_L3V->Draw();
327 if (h_clusterSNR_L456U) h_clusterSNR_L456U->Draw();
329 if (h_clusterSNR_L456V) h_clusterSNR_L456V->Draw();
330
331 float MPVClusterSNRL3U = nan;
332 if (h_clusterSNR_L3U)
333 if (h_clusterSNR_L3U->GetEntries() != 0)
334 MPVClusterSNRL3U = xForMaxY(h_clusterSNR_L3U);
335 float MPVClusterSNRL3V = nan;
336 if (h_clusterSNR_L3V)
337 if (h_clusterSNR_L3V->GetEntries() != 0)
338 MPVClusterSNRL3V = xForMaxY(h_clusterSNR_L3V);
339 float MPVClusterSNRL456U = nan;
340 if (h_clusterSNR_L456U)
341 if (h_clusterSNR_L456U->GetEntries() != 0)
342 MPVClusterSNRL456U = xForMaxY(h_clusterSNR_L456U);
343 float MPVClusterSNRL456V = nan;
344 if (h_clusterSNR_L456V)
345 if (h_clusterSNR_L456V->GetEntries() != 0)
346 MPVClusterSNRL456V = xForMaxY(h_clusterSNR_L456V);
347
348 if (h_clusterSNR_L3U == NULL || h_clusterSNR_L456U == NULL) {
349 B2INFO("Histograms needed for MPV cluster SNR on U side are not found");
350 } else {
351 m_monObj->setVariable("MPVClusterSNRL3U", MPVClusterSNRL3U);
352 m_monObj->setVariable("MPVClusterSNRL456U", MPVClusterSNRL456U);
353 }
354
355 if (h_clusterSNR_L3V == NULL || h_clusterSNR_L456V == NULL) {
356 B2INFO("Histograms needed for MPV cluster SNR on V side are not found");
357 } else {
358 m_monObj->setVariable("MPVClusterSNRL3V", MPVClusterSNRL3V);
359 m_monObj->setVariable("MPVClusterSNRL456V", MPVClusterSNRL456V);
360 }
361
362
363 // MPV SVD cluster time for the clusters on track
364 TH1F* h_clusterTime_L3U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeU3");
365 TH1F* h_clusterTime_L3V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeV3");
366 TH1F* h_clusterTime_L456U = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeU456");
367 TH1F* h_clusterTime_L456V = (TH1F*)findHist("SVDClsTrk/SVDTRK_ClusterTimeV456");
368 TH1F* h_MeanSVD3EventT0 = (TH1F*)findHist("SVDHitTime/SVD3EventT0");
369 TH1F* h_MeanSVD6EventT0 = (TH1F*)findHist("SVDHitTime/SVD6EventT0");
370 TH1F* h_MeanSVDEventT0 = 0x0;
371
372 if (h_MeanSVD3EventT0)
373 h_MeanSVDEventT0 = (TH1F*)h_MeanSVD3EventT0->Clone();
374
376 m_c_MPVTimeClusterOnTrack->Divide(2, 2);
378 if (h_clusterTime_L3U) h_clusterTime_L3U->Draw();
380 if (h_clusterTime_L3V) h_clusterTime_L3V->Draw();
382 if (h_clusterTime_L456U) h_clusterTime_L456U->Draw();
384 if (h_clusterTime_L456V) h_clusterTime_L456V->Draw();
385
386 m_c_MeanSVDEventT0->Clear();
387 m_c_MeanSVDEventT0->Divide(2, 2);
388 m_c_MeanSVDEventT0->cd(1);
389 if (h_MeanSVD3EventT0) h_MeanSVD3EventT0->Draw();
390 m_c_MeanSVDEventT0->cd(2);
391 if (h_MeanSVD6EventT0) h_MeanSVD6EventT0->Draw();
392 m_c_MeanSVDEventT0->cd(3);
393 if (h_MeanSVDEventT0) {
394 if (h_MeanSVD6EventT0)
395 h_MeanSVDEventT0->Add(h_MeanSVD6EventT0);
396 h_MeanSVDEventT0->Draw();
397 }
398
399 float MPVClusterTimeL3U = nan;
400 if (h_clusterTime_L3U)
401 if (h_clusterTime_L3U->GetEntries() != 0)
402 MPVClusterTimeL3U = xForMaxY(h_clusterTime_L3U);
403 float MPVClusterTimeL3V = nan;
404 if (h_clusterTime_L3V)
405 if (h_clusterTime_L3V->GetEntries() != 0)
406 MPVClusterTimeL3V = xForMaxY(h_clusterTime_L3V);
407 float MPVClusterTimeL456U = nan;
408 if (h_clusterTime_L456U)
409 if (h_clusterTime_L456U->GetEntries() != 0)
410 MPVClusterTimeL456U = xForMaxY(h_clusterTime_L456U);
411 float MPVClusterTimeL456V = nan;
412 if (h_clusterTime_L456V)
413 if (h_clusterTime_L456V->GetEntries() != 0)
414 MPVClusterTimeL456V = xForMaxY(h_clusterTime_L456V);
415 float FWHMClusterTimeL3U = nan;
416 if (h_clusterTime_L3U)
417 if (h_clusterTime_L3U->GetEntries() != 0)
418 FWHMClusterTimeL3U = histFWHM(h_clusterTime_L3U);
419 float FWHMClusterTimeL3V = nan;
420 if (h_clusterTime_L3V)
421 if (h_clusterTime_L3V->GetEntries() != 0)
422 FWHMClusterTimeL3V = histFWHM(h_clusterTime_L3V);
423 float FWHMClusterTimeL456U = nan;
424 if (h_clusterTime_L456U)
425 if (h_clusterTime_L456U->GetEntries() != 0)
426 FWHMClusterTimeL456U = histFWHM(h_clusterTime_L456U);
427 float FWHMClusterTimeL456V = nan;
428 if (h_clusterTime_L456V)
429 if (h_clusterTime_L456V->GetEntries() != 0)
430 FWHMClusterTimeL456V = histFWHM(h_clusterTime_L456V);
431
432 float MeanSVD3EventT0 = nan;
433 if (h_MeanSVD3EventT0)
434 if (h_MeanSVD3EventT0->GetEntries() != 0)
435 MeanSVD3EventT0 = xForMaxY(h_MeanSVD3EventT0);
436
437 float MeanSVD6EventT0 = nan;
438 if (h_MeanSVD6EventT0)
439 if (h_MeanSVD6EventT0->GetEntries() != 0)
440 MeanSVD6EventT0 = xForMaxY(h_MeanSVD6EventT0);
441
442 float MeanSVDEventT0 = nan;
443 if (h_MeanSVDEventT0)
444 if (h_MeanSVDEventT0->GetEntries() != 0)
445 MeanSVDEventT0 = xForMaxY(h_MeanSVDEventT0);
446
447 if (h_clusterTime_L3U == NULL || h_clusterTime_L456U == NULL) {
448 B2INFO("Histograms needed for MPV cluster time on U side are not found");
449 } else {
450 m_monObj->setVariable("MPVClusterTimeL3U", MPVClusterTimeL3U);
451 m_monObj->setVariable("MPVClusterTimeL456U", MPVClusterTimeL456U);
452 m_monObj->setVariable("FWHMClusterTimeL3U", FWHMClusterTimeL3U);
453 m_monObj->setVariable("FWHMClusterTimeL456U", FWHMClusterTimeL456U);
454 }
455
456 if (h_clusterTime_L3V == NULL || h_clusterTime_L456V == NULL) {
457 B2INFO("Histograms needed for MPV cluster time on V side are not found");
458 } else {
459 m_monObj->setVariable("MPVClusterTimeL3V", MPVClusterTimeL3V);
460 m_monObj->setVariable("MPVClusterTimeL456V", MPVClusterTimeL456V);
461 m_monObj->setVariable("FWHMClusterTimeL3V", FWHMClusterTimeL3V);
462 m_monObj->setVariable("FWHMClusterTimeL456V", FWHMClusterTimeL456V);
463 }
464
465 if (h_MeanSVD3EventT0 == NULL) {
466 B2INFO("Histograms needed for SVD Event T0 (3 samples) not found");
467 } else {
468 m_monObj->setVariable("MeanSVD3EventT0", MeanSVD3EventT0);
469 }
470
471 if (h_MeanSVD6EventT0 == NULL) {
472 B2INFO("Histograms needed for SVD Event T0 (6 samples) not found");
473 } else {
474 m_monObj->setVariable("MeanSVD6EventT0", MeanSVD6EventT0);
475 }
476
477 if (h_MeanSVDEventT0 == NULL) {
478 B2INFO("Histograms needed for SVD Event T0 (all samples) not found");
479 } else {
480 m_monObj->setVariable("MeanSVDEventT0", MeanSVDEventT0);
481 }
482
483 // average maxBin for clusters on track
484 TH1F* h_maxBinU = (TH1F*)findHist("SVDClsTrk/SVDTRK_StripMaxBinUAll");
485 TH1F* h_maxBinV = (TH1F*)findHist("SVDClsTrk/SVDTRK_StripMaxBinVAll");
486
488 m_c_avgMaxBinClusterOnTrack->Divide(2, 1);
490 if (h_maxBinU) h_maxBinU->Draw();
492 if (h_maxBinV) h_maxBinV->Draw();
493
494 if (h_maxBinU == NULL) {
495 B2INFO("Histogram needed for Average MaxBin on U side is not found");
496 } else {
497 float avgMaxBinU = h_maxBinU->GetMean();
498 m_monObj->setVariable("avgMaxBinU", avgMaxBinU);
499 }
500
501 if (h_maxBinV == NULL) {
502 B2INFO("Histogram needed for Average MaxBin on V side is not found");
503 } else {
504 float avgMaxBinV = h_maxBinV->GetMean();
505 m_monObj->setVariable("avgMaxBinV", avgMaxBinV);
506 }
507
508 // Cluster on track ladder
509 for (const auto& it : ladderLabel) {
510 string sensorDescr = it;
511 int layer = 0;
512 int sensor = 0;
513 sscanf(it.c_str(), "L%d.X.%d", &layer, &sensor);
514
515 TString name = Form("SVDClsTrk/SVDTRK_ClusterCharge_L%d.x.%d", layer, sensor);
516 TString title = Form("MPVClusterCharge_L%d.x.%d", layer, sensor);
517 TH1F* h_clusterCharge = (TH1F*)findHist(name.Data());
518 float MPVClusterCharge = nan;
519 if (h_clusterCharge)
520 if (h_clusterCharge->GetEntries() != 0)
521 MPVClusterCharge = xForMaxY(h_clusterCharge);
522
523 if (h_clusterCharge == NULL) {
524 B2INFO("Histograms needed for cluster charge not found");
525 } else {
526 m_monObj->setVariable(title.Data(), MPVClusterCharge);
527 }
528
529 name = Form("SVDClsTrk/SVDTRK_ClusterSNR_L%d.x.%d", layer, sensor);
530 title = Form("MPVClusterSNR_L%d.x.%d", layer, sensor);
531 TH1F* h_clusterSNR = (TH1F*)findHist(name.Data());
532 float MPVClusterSNR = nan;
533 if (h_clusterSNR)
534 if (h_clusterSNR->GetEntries() != 0)
535 MPVClusterSNR = xForMaxY(h_clusterSNR);
536
537 if (h_clusterSNR == NULL) {
538 B2INFO("Histograms needed for cluster SNR not found");
539 } else {
540 m_monObj->setVariable(title.Data(), MPVClusterSNR);
541 }
542 }
543
544 // Cluster on track peculiar sensors
545 for (const auto& it : m_listOfSensorsToMonitor) {
546 string sensorDescr = it;
547 string valueLabel = it;
548 replace(sensorDescr.begin(), sensorDescr.end(), '.', '_');
549 valueLabel.erase(remove(valueLabel.begin(), valueLabel.end(), '.'), valueLabel.end());
550
551 TString name = Form("SVDClsTrk/SVDTRK_%s_ClusterChargeU", sensorDescr.c_str());
552 TString title = Form("MPVClusterChargeL%sU", valueLabel.c_str());
553 TString title1 = "";
554 TH1F* h_clusterCharge = (TH1F*)findHist(name.Data());
555 float MPVClusterCharge = nan;
556 if (h_clusterCharge)
557 if (h_clusterCharge->GetEntries() != 0)
558 MPVClusterCharge = xForMaxY(h_clusterCharge);
559
560 if (h_clusterCharge == NULL) {
561 B2INFO("Histograms needed for clusterU charge not found");
562 } else {
563 m_monObj->setVariable(title.Data(), MPVClusterCharge);
564 }
565
566 name = Form("SVDClsTrk/SVDTRK_%s_ClusterChargeV", sensorDescr.c_str());
567 title = Form("MPVClusterChargeL%sV", valueLabel.c_str());
568 h_clusterCharge = (TH1F*)findHist(name.Data());
569 MPVClusterCharge = nan;
570 if (h_clusterCharge)
571 if (h_clusterCharge->GetEntries() != 0)
572 MPVClusterCharge = xForMaxY(h_clusterCharge);
573
574 if (h_clusterCharge == NULL) {
575 B2INFO("Histograms needed for clusterV charge not found");
576 } else {
577 m_monObj->setVariable(title.Data(), MPVClusterCharge);
578 }
579
580 name = Form("SVDClsTrk/SVDTRK_%s_ClusterSNRU", sensorDescr.c_str());
581 title = Form("MPVClusterSNRL%sU", valueLabel.c_str());
582 TH1F* h_clusterSNR = (TH1F*)findHist(name.Data());
583 float MPVClusterSNR = nan;
584 if (h_clusterSNR)
585 if (h_clusterSNR->GetEntries() != 0)
586 MPVClusterSNR = xForMaxY(h_clusterSNR);
587
588 if (h_clusterSNR == NULL) {
589 B2INFO("Histograms needed for clusterU SNR not found");
590 } else {
591 m_monObj->setVariable(title.Data(), MPVClusterSNR);
592 }
593
594 name = Form("SVDClsTrk/SVDTRK_%s_ClusterSNRV", sensorDescr.c_str());
595 title = Form("MPVClusterSNRL%sV", valueLabel.c_str());
596 h_clusterSNR = (TH1F*)findHist(name.Data());
597 MPVClusterSNR = nan;
598 if (h_clusterSNR)
599 if (h_clusterSNR->GetEntries() != 0)
600 MPVClusterSNR = xForMaxY(h_clusterSNR);
601
602 if (h_clusterSNR == NULL) {
603 B2INFO("Histograms needed for clusterV SNR not found");
604 } else {
605 m_monObj->setVariable(title.Data(), MPVClusterSNR);
606 }
607
608 name = Form("SVDClsTrk/SVDTRK_%s_ClusterTimeU", sensorDescr.c_str());
609 title = Form("MPVClusterTimeL%sU", valueLabel.c_str());
610 title1 = Form("FWHMClusterTimeL%sU", valueLabel.c_str());
611 TH1F* h_clusterTime = (TH1F*)findHist(name.Data());
612 float MPVClusterTime = nan;
613 float FWHMClusterTime = nan;
614 if (h_clusterTime)
615 if (h_clusterTime->GetEntries() != 0) {
616 MPVClusterTime = xForMaxY(h_clusterTime);
617 FWHMClusterTime = histFWHM(h_clusterTime);
618 }
619
620 if (h_clusterTime == NULL) {
621 B2INFO("Histograms needed for clusterU time not found");
622 } else {
623 m_monObj->setVariable(title.Data(), MPVClusterTime);
624 m_monObj->setVariable(title1.Data(), FWHMClusterTime);
625 }
626
627 name = Form("SVDClsTrk/SVDTRK_%s_ClusterTimeV", sensorDescr.c_str());
628 title = Form("MPVClusterTimeL%sV", valueLabel.c_str());
629 title1 = Form("FWHMClusterTimeL%sV", valueLabel.c_str());
630 h_clusterTime = (TH1F*)findHist(name.Data());
631 MPVClusterTime = nan;
632 FWHMClusterTime = nan;
633 if (h_clusterTime)
634 if (h_clusterTime->GetEntries() != 0) {
635 MPVClusterTime = xForMaxY(h_clusterTime);
636 FWHMClusterTime = histFWHM(h_clusterTime);
637 }
638
639 if (h_clusterTime == NULL) {
640 B2INFO("Histograms needed for clusterU time not found");
641 } else {
642 m_monObj->setVariable(title.Data(), MPVClusterTime);
643 m_monObj->setVariable(title1.Data(), FWHMClusterTime);
644 }
645 }
646
647 B2INFO("DQMHistAnalysisSVDGeneral: endRun called");
648}
649
650
652{
653 B2INFO("DQMHistAnalysisSVDOnMiraBelle: terminate called");
654}
655
656std::pair<float, float> DQMHistAnalysisSVDOnMiraBelleModule::avgOccupancyUV(TH1F* hU, TH1F* hV, int nEvents,
657 int layer, int ladder, int sensor) const
658{
659 int nStripsV = -1;
660 if (layer == 3) {
661 nStripsV = 768;
662 } else if (layer >= 4 && layer <= 6) {
663 nStripsV = 512;
664 } else {
665 B2DEBUG(20, "Layer out of range [3,6].");
666 }
667 std::pair<float, float> avgOffOccUV(0.0, 0.0);
668
669 int minLayer = (layer != -1) ? layer : m_gTools->getFirstSVDLayer();
670 int maxLayer = (layer != -1) ? layer : m_gTools->getLastSVDLayer();
671 int sensorsN = 0;
672
673 if (ladder == 0) ladder = -1;
674
675 for (int layerId = minLayer; layerId < maxLayer + 1; ++layerId) {
676 int minLadder = (ladder != -1) ? ladder : 1;
677 int maxLadder = (ladder != -1) ? ladder : getNumberOfLadders(layerId);
678
679 int minSensor = (sensor != -1) ? sensor : 1;
680 int maxSensor = (sensor != -1) ? sensor : getNumberOfSensors(layerId);
681
682 for (int sensorId = minSensor; sensorId < maxSensor + 1; ++sensorId) {
683
684 for (int ladderId = minLadder; ladderId < maxLadder + 1; ++ladderId) {
685 int bin = m_gTools->getSVDSensorIndex(layerId, ladderId, sensorId) + 1;
686
687 avgOffOccUV.first += hU->GetBinContent(bin) / 768 * 100;
688 avgOffOccUV.second += hV->GetBinContent(bin) / nStripsV * 100;
689 sensorsN++;
690 }
691 }
692 }
693
694 avgOffOccUV.first /= (sensorsN * nEvents);
695 avgOffOccUV.second /= (sensorsN * nEvents);
696
697 return avgOffOccUV;
698}
699
700std::pair<float, float> DQMHistAnalysisSVDOnMiraBelleModule::avgOccupancyGrpId0UV(int iLayer, int nEvents) const
701{
702 int nStripsV = -1;
703 if (iLayer == 3) {
704 nStripsV = 768;
705 } else if (iLayer >= 4 && iLayer <= 6) {
706 nStripsV = 512;
707 } else {
708 B2DEBUG(20, "Layer out of range [3,6].");
709 }
710
711 Int_t nStripsU = 768;
712
713 std::vector<float> avgOffOccU;
714 std::vector<float> avgOffOccV;
715
716 for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
717 int tmp_layer = m_SVDModules[i].getLayerNumber();
718 int tmp_ladder = m_SVDModules[i].getLadderNumber();
719 int tmp_sensor = m_SVDModules[i].getSensorNumber();
720
721 TString tmpnameGrpId0U = Form("SVDExpReco/SVDDQM_%d_%d_%d_StripCountGroupId0U", tmp_layer, tmp_ladder, tmp_sensor);
722 TH1F* htmpU = (TH1F*)findHist(tmpnameGrpId0U.Data());
723 if (htmpU == NULL) {
724 B2INFO("Occupancy U histogram for group Id0 not found");
725 } else {
726 if (tmp_layer == iLayer)
727 avgOffOccU.push_back(htmpU->GetEntries() / nStripsU / nEvents * 100);
728 }
729
730 TString tmpnameGrpId0V = Form("SVDExpReco/SVDDQM_%d_%d_%d_StripCountGroupId0V", tmp_layer, tmp_ladder, tmp_sensor);
731 TH1F* htmpV = (TH1F*)findHist(tmpnameGrpId0V.Data());
732 if (htmpV == NULL) {
733 B2INFO("Occupancy V histogram for group Id0 not found");
734 } else {
735 if (tmp_layer == iLayer)
736 avgOffOccV.push_back(htmpV->GetEntries() / nStripsV / nEvents * 100);
737 }
738 }
739
740 std::pair<float, float> avgOffOccUV(0., 0.);
741
742 avgOffOccUV.first = accumulate(avgOffOccU.begin(), avgOffOccU.end(), 0.0);
743 avgOffOccUV.first /= float(avgOffOccU.size());
744
745 avgOffOccUV.second = accumulate(avgOffOccV.begin(), avgOffOccV.end(), 0.0);
746 avgOffOccUV.second /= float(avgOffOccV.size());
747
748 return avgOffOccUV;
749}
750
751std::pair<float, float> DQMHistAnalysisSVDOnMiraBelleModule::avgEfficiencyUV(TH2F* hMCU, TH2F* hMCV, TH2F* hFTU, TH2F* hFTV,
752 int layer,
753 int ladder, int sensor) const
754{
755 float nan = numeric_limits<float>::quiet_NaN();
756 std::pair<float, float> avgEffUV(0.0, 0.0);
757 std::pair<float, float> sumMatchedClustersUV(0.0, 0.0);
758 std::pair<float, float> sumFoundTracksUV(0.0, 0.0);
759
760 int minLayer = (layer != -1) ? layer : m_gTools->getFirstSVDLayer();
761 int maxLayer = (layer != -1) ? layer : m_gTools->getLastSVDLayer();
762
763 if (ladder == 0) ladder = -1;
764
765 for (int layerId = minLayer; layerId < maxLayer + 1; ++layerId) {
766 int minLadder = (ladder != -1) ? ladder : 1;
767 int maxLadder = (ladder != -1) ? ladder : getNumberOfLadders(layerId);
768
769 int minSensor = (sensor != -1) ? sensor : 1;
770 int maxSensor = (sensor != -1) ? sensor : getNumberOfSensors(layerId);
771
772 for (int sensorId = minSensor; sensorId < maxSensor + 1; ++sensorId) {
773
774 for (int ladderId = minLadder; ladderId < maxLadder + 1; ++ladderId) {
775 int binY = findBinY(layerId, sensorId);
776 int binXY = hMCV->FindBin(ladderId, binY);
777
778 sumMatchedClustersUV.first += hMCU->GetBinContent(binXY);
779 sumMatchedClustersUV.second += hMCV->GetBinContent(binXY);
780 sumFoundTracksUV.first += hFTU->GetBinContent(binXY);
781 sumFoundTracksUV.second += hFTV->GetBinContent(binXY);
782 }
783 }
784 }
785
786 if (sumFoundTracksUV.first > 0) {
787 avgEffUV.first = sumMatchedClustersUV.first / sumFoundTracksUV.first * 100;
788 } else {
789 avgEffUV.first = nan;
790 }
791
792 if (sumFoundTracksUV.second > 0) {
793 avgEffUV.second = sumMatchedClustersUV.second / sumFoundTracksUV.second * 100;
794 } else {
795 avgEffUV.second = nan;
796 }
797
798 return avgEffUV;
799}
800
801
802void DQMHistAnalysisSVDOnMiraBelleModule::addVariable(string name, pair<float, float>& varUV)
803{
804 size_t pos = name.find("UV");
805
806 if (pos != string::npos)
807 name.replace(pos, 2, "");
808
809 m_monObj->setVariable(Form("%sU", name.c_str()), varUV.first);
810 m_monObj->setVariable(Form("%sV", name.c_str()), varUV.second);
811}
812
814{
815 int maxY = h->GetMaximumBin();
816 float xMaxY = h->GetXaxis()->GetBinCenter(maxY);
817 return xMaxY;
818}
819
821{
822 int bin1 = h->FindFirstBinAbove(h->GetMaximum() / 2);
823 int bin2 = h->FindLastBinAbove(h->GetMaximum() / 2);
824 float fwhm = h->GetBinCenter(bin2) - h->GetBinCenter(bin1);
825 return fwhm;
826}
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
DQMHistAnalysisModule()
Constructor / Destructor.
void initialize() override final
Module function initialize.
float histFWHM(TH1F *h) const
Calculate full width at half maximum of histogram.
TCanvas * m_c_MPVTimeClusterOnTrack
time for clusters on track
std::vector< VxdID > m_SVDModules
IDs of all SVD Modules to iterate over.
const VXD::GeoTools * m_gTools
geometrical tool pointer
void addVariable(std::string name, std::pair< float, float > &varUV)
Add variable to object monitoring.
TCanvas * m_c_MPVChargeClusterOnTrack
charge for clusters on track
DBObjPtr< SVDDQMPlotsConfiguration > m_svdPlotsConfig
SVD DQM plots configuration.
MonitoringObject * m_monObj
Monitoring Object to be produced by this module, which contain defined canvases and monitoring variab...
TCanvas * m_c_avgEfficiency
List of canvases to be added to MonitoringObject.
void terminate() override final
Module function terminate.
Int_t findBinY(Int_t layer, Int_t sensor) const
find the Y bin given the layer and sensor number
void event() override final
Module function event.
std::pair< float, float > avgOccupancyGrpId0UV(int iLayer, int nEvents) const
Calculate avg offline occupancy for specified layer for time group id = 0.
std::pair< float, float > avgOccupancyUV(TH1F *hU, TH1F *hV, int nEvents, int layer=-1, int ladder=-1, int sensor=-1) const
Calculate avg offline occupancy for one specific sensor, especially.
std::pair< float, float > avgEfficiencyUV(TH2F *hMCU, TH2F *hMCV, TH2F *hFTU, TH2F *hFTV, int layer=-1, int ladder=-1, int sensor=-1) const
Calculate avg efficiency for specified sensors.
TCanvas * m_c_MPVSNRClusterOnTrack
SNR for clusters on track.
Int_t getNumberOfLadders(Int_t layer) const
get number of ladders per layer
void endRun() override final
Module function endRun.
float xForMaxY(TH1F *h) const
Calculate abscissa of max Y bin.
TCanvas * m_c_avgOffOccupancy
number of ZS5 fired strips
void beginRun() override final
Module function beginRun.
TCanvas * m_c_avgMaxBinClusterOnTrack
average number of the APV sample which corresponds to the maximum amplitude for clusters on track
std::vector< std::string > m_listOfSensorsToMonitor
list of sensor to monitor (Charge, SNR, time; U/V) taken from DB (payload)
Int_t getNumberOfSensors(Int_t layer) const
get number of sensors per layer
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition Module.h:80
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
Definition GeoCache.h:38
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a reference 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
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition GeoCache.h:141
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:32
#define SetVariable(x)
set variable to mirabelle for a given member
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.
STL namespace.