Belle II Software  release-05-01-25
SVDDQMExpressRecoModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Peter Kodys, Giulia Casarosa, Giuliana Rizzo *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include "svd/modules/svdDQM/SVDDQMExpressRecoModule.h"
12 
13 #include <hlt/softwaretrigger/core/FinalTriggerDecisionCalculator.h>
14 #include <framework/datastore/StoreObjPtr.h>
15 #include <framework/datastore/StoreArray.h>
16 #include <framework/dataobjects/EventMetaData.h>
17 
18 #include <svd/dataobjects/SVDShaperDigit.h>
19 #include <svd/dataobjects/SVDCluster.h>
20 
21 #include <vxd/geometry/SensorInfoBase.h>
22 #include <vxd/geometry/GeoTools.h>
23 
24 #include <boost/format.hpp>
25 
26 #include "TDirectory.h"
27 
28 using namespace std;
29 using boost::format;
30 using namespace Belle2;
31 using namespace SoftwareTrigger;
32 
33 //-----------------------------------------------------------------
34 // Register the Module
35 //-----------------------------------------------------------------
36 REG_MODULE(SVDDQMExpressReco)
37 
38 
39 //-----------------------------------------------------------------
40 // Implementation
41 //-----------------------------------------------------------------
42 
44 {
45  //Set module properties
46  setDescription("SVD DQM module for Express Reco"
47  "Recommended Number of events for monitor is 40 kEvents or more to fill all histograms "
48  "Container for histograms for off-line analysis with any granulation base on request "
49  "Call all histograms set ShowAllHistos=1 ."
50  );
51 
52  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
53  addParam("offlineZSShaperDigits", m_storeSVDShaperDigitsName, "ShaperDigits StoreArray name", std::string("SVDShaperDigitsZS5"));
54  addParam("ShaperDigits", m_storeNoZSSVDShaperDigitsName, "not zero-suppressed ShaperDigits StoreArray name",
55  std::string("SVDShaperDigits"));
56  addParam("skipHLTRejectedEvents", m_skipRejectedEvents, "If TRUE skip events rejected by HLT", bool(true));
57  addParam("ShowAllHistos", m_ShowAllHistos, "Flag to show all histos in DQM, default = 0 ", int(0));
58  addParam("desynchronizeSVDTime", m_desynchSVDTime,
59  "if TRUE (default is FALSE): svdTime back in SVD time reference", bool(false));
60  addParam("CutSVDCharge", m_CutSVDCharge,
61  "cut for accepting to hitmap histogram, using strips only", float(0));
62  addParam("CutSVDClusterCharge", m_CutSVDClusterCharge,
63  "cut for accepting clusters to hitmap histogram", float(0));
64  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
65  std::string("SVDExpReco"));
66 
67  m_histoList = new TList();
68 }
69 
70 
71 SVDDQMExpressRecoModule::~SVDDQMExpressRecoModule()
72 {
73 }
74 
75 //------------------------------------------------------------------
76 // Function to define histograms
77 //-----------------------------------------------------------------
78 
79 void SVDDQMExpressRecoModule::defineHisto()
80 {
81  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
82  if (gTools->getNumberOfLayers() == 0) {
83  B2FATAL("Missing geometry for VXD, check steering file.");
84  }
85  if (gTools->getNumberOfSVDLayers() == 0) {
86  B2WARNING("Missing geometry for SVD, SVD-DQM is skiped.");
87  return;
88  }
89 
90  // Create a separate histogram directories and cd into it.
91  TDirectory* oldDir = gDirectory;
92  if (m_histogramDirectoryName != "") {
93  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
94  oldDir->cd(m_histogramDirectoryName.c_str());
95  }
96 
97  // basic constants presets:
98  int nSVDSensors = gTools->getNumberOfSVDSensors();
99  int nSVDChips = gTools->getTotalSVDChips();
100 
101  // number of events counter
102  m_nEvents = new TH1F("SVDDQM_nEvents", "SVD Number of Events", 1, -0.5, 0.5);
103  m_nEvents->GetYaxis()->SetTitle("N events");
104  m_histoList->Add(m_nEvents);
105 
106  // Create basic histograms:
107  // basic counters per sensor:
108  m_hitMapCountsU = new TH1F("SVDDQM_StripCountsU", "SVD Integrated Number of ZS5 Fired U-Strips per sensor",
109  nSVDSensors, 0, nSVDSensors);
110  m_hitMapCountsU->GetXaxis()->SetTitle("Sensor ID");
111  m_hitMapCountsU->GetYaxis()->SetTitle("counts");
112  m_histoList->Add(m_hitMapCountsU);
113  m_hitMapCountsV = new TH1F("SVDDQM_StripCountsV", "SVD Integrated Number of ZS5 Fired V-Strips per sensor",
114  nSVDSensors, 0, nSVDSensors);
115  m_hitMapCountsV->GetXaxis()->SetTitle("Sensor ID");
116  m_hitMapCountsV->GetYaxis()->SetTitle("counts");
117  m_histoList->Add(m_hitMapCountsV);
118  m_hitMapClCountsU = new TH1F("SVDDQM_ClusterCountsU", "SVD Integrated Number of U-Clusters per sensor",
119  nSVDSensors, 0, nSVDSensors);
120  m_hitMapClCountsU->GetXaxis()->SetTitle("Sensor ID");
121  m_hitMapClCountsU->GetYaxis()->SetTitle("counts");
122  m_histoList->Add(m_hitMapClCountsU);
123  m_hitMapClCountsV = new TH1F("SVDDQM_ClusterCountsV", "SVD Integrated Number of V-Clusters per sensor",
124  nSVDSensors, 0, nSVDSensors);
125  m_hitMapClCountsV->GetXaxis()->SetTitle("Sensor ID");
126  m_hitMapClCountsV->GetYaxis()->SetTitle("counts");
127  m_histoList->Add(m_hitMapClCountsV);
128  for (int i = 0; i < nSVDSensors; i++) {
129  VxdID id = gTools->getSensorIDFromSVDIndex(i);
130  int iLayer = id.getLayerNumber();
131  int iLadder = id.getLadderNumber();
132  int iSensor = id.getSensorNumber();
133  TString AxisTicks = Form("%i_%i_%i", iLayer, iLadder, iSensor);
134  m_hitMapCountsU->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
135  m_hitMapCountsV->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
136  m_hitMapClCountsU->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
137  m_hitMapClCountsV->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
138  }
139 
140  // basic counters per chip:
141  m_hitMapCountsChip = new TH1F("SVDDQM_StripCountsChip", "SVD Integrated Number of ZS5 Fired Strips per chip",
142  nSVDChips, 0, nSVDChips);
143  m_hitMapCountsChip->GetXaxis()->SetTitle("Chip ID");
144  m_hitMapCountsChip->GetYaxis()->SetTitle("counts");
145  m_histoList->Add(m_hitMapCountsChip);
146  m_hitMapClCountsChip = new TH1F("SVDDQM_ClusterCountsChip", "SVD Integrated Number of Clusters per chip",
147  nSVDChips, 0, nSVDChips);
148  m_hitMapClCountsChip->GetXaxis()->SetTitle("Chip ID");
149  m_hitMapClCountsChip->GetYaxis()->SetTitle("counts");
150  m_histoList->Add(m_hitMapClCountsChip);
151 
152  m_firedU = new TH1F*[nSVDSensors];
153  m_firedV = new TH1F*[nSVDSensors];
154  m_clustersU = new TH1F*[nSVDSensors];
155  m_clustersV = new TH1F*[nSVDSensors];
156 
157 
158  m_clusterChargeU = new TH1F*[nSVDSensors];
159  m_clusterChargeV = new TH1F*[nSVDSensors];
160  m_clusterSNRU = new TH1F*[nSVDSensors];
161  m_clusterSNRV = new TH1F*[nSVDSensors];
162  m_stripSignalU = new TH1F*[nSVDSensors];
163  m_stripSignalV = new TH1F*[nSVDSensors];
164  m_stripCountU = new TH1F*[nSVDSensors];
165  m_stripCountV = new TH1F*[nSVDSensors];
166  m_onlineZSstripCountU = new TH1F*[nSVDSensors];
167  m_onlineZSstripCountV = new TH1F*[nSVDSensors];
168  m_clusterSizeU = new TH1F*[nSVDSensors];
169  m_clusterSizeV = new TH1F*[nSVDSensors];
170  m_clusterTimeU = new TH1F*[nSVDSensors];
171  m_clusterTimeV = new TH1F*[nSVDSensors];
172 
173  int ChargeBins = 80;
174  float ChargeMax = 80;
175  int SNRBins = 50;
176  float SNRMax = 100;
177  int TimeBins = 300;
178  float TimeMin = -150;
179  float TimeMax = 150;
180 
181  int MaxBinBins = 6;
182  int MaxBinMax = 6;
183 
184  TString refFrame = "in FTSW reference";
185  if (m_desynchSVDTime)
186  refFrame = "in SVD reference";
187 
188 
189  //----------------------------------------------------------------
190  // Charge of clusters for all sensors
191  //----------------------------------------------------------------
192  string name = str(format("SVDDQM_ClusterChargeUAll"));
193  string title = str(format("SVD U-Cluster Charge for all sensors"));
194  m_clusterChargeUAll = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
195  m_clusterChargeUAll->GetXaxis()->SetTitle("cluster charge [ke-]");
196  m_clusterChargeUAll->GetYaxis()->SetTitle("count");
197  m_histoList->Add(m_clusterChargeUAll);
198  name = str(format("SVDDQM_ClusterChargeVAll"));
199  title = str(format("SVD V-Cluster Charge for all sensors"));
200  m_clusterChargeVAll = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
201  m_clusterChargeVAll->GetXaxis()->SetTitle("cluster charge [ke-]");
202  m_clusterChargeVAll->GetYaxis()->SetTitle("count");
203  m_histoList->Add(m_clusterChargeVAll);
204  //----------------------------------------------------------------
205  // Charge of clusters for L3/L456 sensors
206  //----------------------------------------------------------------
207  name = str(format("SVDDQM_ClusterChargeU3"));
208  title = str(format("SVD U-Cluster Charge for layer 3 sensors"));
209  m_clusterChargeU3 = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
210  m_clusterChargeU3->GetXaxis()->SetTitle("cluster charge [ke-]");
211  m_clusterChargeU3->GetYaxis()->SetTitle("count");
212  m_histoList->Add(m_clusterChargeU3);
213  name = str(format("SVDDQM_ClusterChargeV3"));
214  title = str(format("SVD V-Cluster Charge for layer 3 sensors"));
215  m_clusterChargeV3 = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
216  m_clusterChargeV3->GetXaxis()->SetTitle("cluster charge [ke-]");
217  m_clusterChargeV3->GetYaxis()->SetTitle("count");
218  m_histoList->Add(m_clusterChargeV3);
219 
220  name = str(format("SVDDQM_ClusterChargeU456"));
221  title = str(format("SVD U-Cluster Charge for layers 4,5,6 sensors"));
222  m_clusterChargeU456 = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
223  m_clusterChargeU456->GetXaxis()->SetTitle("cluster charge [ke-]");
224  m_clusterChargeU456->GetYaxis()->SetTitle("count");
225  m_histoList->Add(m_clusterChargeU456);
226 
227  name = str(format("SVDDQM_ClusterChargeV456"));
228  title = str(format("SVD V-Cluster Charge for layers 4,5,6 sensors"));
229  m_clusterChargeV456 = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
230  m_clusterChargeV456->GetXaxis()->SetTitle("cluster charge [ke-]");
231  m_clusterChargeV456->GetYaxis()->SetTitle("count");
232  m_histoList->Add(m_clusterChargeV456);
233 
234  //----------------------------------------------------------------
235  // SNR of clusters for all sensors
236  //----------------------------------------------------------------
237  name = str(format("SVDDQM_ClusterSNRUAll"));
238  title = str(format("SVD U-Cluster SNR for all sensors"));
239  m_clusterSNRUAll = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax); // max = ~ 60
240  m_clusterSNRUAll->GetXaxis()->SetTitle("cluster SNR");
241  m_clusterSNRUAll->GetYaxis()->SetTitle("count");
242  m_histoList->Add(m_clusterSNRUAll);
243  name = str(format("SVDDQM_ClusterSNRVAll"));
244  title = str(format("SVD V-Cluster SNR for all sensors"));
245  m_clusterSNRVAll = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
246  m_clusterSNRVAll->GetXaxis()->SetTitle("cluster SNR");
247  m_clusterSNRVAll->GetYaxis()->SetTitle("count");
248  m_histoList->Add(m_clusterSNRVAll);
249  //----------------------------------------------------------------
250  // SNR of clusters for L3/L456 sensors
251  //----------------------------------------------------------------
252  name = str(format("SVDDQM_ClusterSNRU3"));
253  title = str(format("SVD U-Cluster SNR for layer 3 sensors"));
254  m_clusterSNRU3 = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
255  m_clusterSNRU3->GetXaxis()->SetTitle("cluster SNR");
256  m_clusterSNRU3->GetYaxis()->SetTitle("count");
257  m_histoList->Add(m_clusterSNRU3);
258  name = str(format("SVDDQM_ClusterSNRV3"));
259  title = str(format("SVD V-Cluster SNR for layer 3 sensors"));
260  m_clusterSNRV3 = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
261  m_clusterSNRV3->GetXaxis()->SetTitle("cluster SNR");
262  m_clusterSNRV3->GetYaxis()->SetTitle("count");
263  m_histoList->Add(m_clusterSNRV3);
264 
265  name = str(format("SVDDQM_ClusterSNRU456"));
266  title = str(format("SVD U-Cluster SNR for layers 4,5,6 sensors"));
267  m_clusterSNRU456 = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
268  m_clusterSNRU456->GetXaxis()->SetTitle("cluster SNR");
269  m_clusterSNRU456->GetYaxis()->SetTitle("count");
270  m_histoList->Add(m_clusterSNRU456);
271  name = str(format("SVDDQM_ClusterSNRV456"));
272  title = str(format("SVD V-Cluster SNR for layers 4,5,6 sensors"));
273  m_clusterSNRV456 = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
274  m_clusterSNRV456->GetXaxis()->SetTitle("cluster SNR");
275  m_clusterSNRV456->GetYaxis()->SetTitle("count");
276  m_histoList->Add(m_clusterSNRV456);
277  //----------------------------------------------------------------
278  // Cluster time distribution for all sensors
279  //----------------------------------------------------------------
280  TString Name = "SVDDQM_ClusterTimeUAll";
281  TString Title = Form("SVD U-Cluster Time %s for all sensors", refFrame.Data());
282  m_clusterTimeUAll = new TH1F(Name.Data(), Title.Data(), TimeBins, TimeMin, TimeMax);
283  m_clusterTimeUAll->GetXaxis()->SetTitle("cluster time (ns)");
284  m_clusterTimeUAll->GetYaxis()->SetTitle("count");
285  m_histoList->Add(m_clusterTimeUAll);
286  Name = "SVDDQM_ClusterTimeVAll";
287  Title = Form("SVD V-Cluster Time %s for all sensors", refFrame.Data());
288  m_clusterTimeVAll = new TH1F(Name.Data(), Title.Data(), TimeBins, TimeMin, TimeMax);
289  m_clusterTimeVAll->GetXaxis()->SetTitle("cluster time (ns)");
290  m_clusterTimeVAll->GetYaxis()->SetTitle("count");
291  m_histoList->Add(m_clusterTimeVAll);
292  //----------------------------------------------------------------
293  // Time of clusters for L3/L456 sensors
294  //----------------------------------------------------------------
295  Name = "SVDDQM_ClusterTimeU3";
296  Title = Form("SVD U-Cluster Time %s for layer 3 sensors", refFrame.Data());
297  m_clusterTimeU3 = new TH1F(Name.Data(), Title.Data(), TimeBins, TimeMin, TimeMax);
298  m_clusterTimeU3->GetXaxis()->SetTitle("cluster time (ns)");
299  m_clusterTimeU3->GetYaxis()->SetTitle("count");
300  m_histoList->Add(m_clusterTimeU3);
301  name = str(format("SVDDQM_ClusterTimeV3"));
302  Title = Form("SVD V-Cluster Time %s for layer 3 sensors", refFrame.Data());
303  m_clusterTimeV3 = new TH1F(name.c_str(), Title.Data(), TimeBins, TimeMin, TimeMax);
304  m_clusterTimeV3->GetXaxis()->SetTitle("cluster time (ns)");
305  m_clusterTimeV3->GetYaxis()->SetTitle("count");
306  m_histoList->Add(m_clusterTimeV3);
307 
308  name = str(format("SVDDQM_ClusterTimeU456"));
309  Title = Form("SVD U-Cluster Time %s for layers 4,5,6 sensors", refFrame.Data());
310  m_clusterTimeU456 = new TH1F(name.c_str(), Title.Data(), TimeBins, TimeMin, TimeMax);
311  m_clusterTimeU456->GetXaxis()->SetTitle("cluster time (ns)");
312  m_clusterTimeU456->GetYaxis()->SetTitle("count");
313  m_histoList->Add(m_clusterTimeU456);
314  name = str(format("SVDDQM_ClusterTimeV456"));
315  Title = Form("SVD V-Cluster Time %s for layers 4,5,6 sensors", refFrame.Data());
316  m_clusterTimeV456 = new TH1F(name.c_str(), Title.Data(), TimeBins, TimeMin, TimeMax);
317  m_clusterTimeV456->GetXaxis()->SetTitle("cluster time (ns)");
318  m_clusterTimeV456->GetYaxis()->SetTitle("count");
319  m_histoList->Add(m_clusterTimeV456);
320 
321  //----------------------------------------------------------------
322  // MaxBin of strips for all sensors (offline ZS)
323  //----------------------------------------------------------------
324  name = str(format("SVDDQM_StripMaxBinUAll"));
325  title = str(format("SVD U-Strip MaxBin for all sensors"));
326  m_stripMaxBinUAll = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
327  m_stripMaxBinUAll->GetXaxis()->SetTitle("max bin");
328  m_stripMaxBinUAll->GetYaxis()->SetTitle("count");
329  m_histoList->Add(m_stripMaxBinUAll);
330  name = str(format("SVDDQM_StripMaxBinVAll"));
331  title = str(format("SVD V-Strip MaxBin for all sensors"));
332  m_stripMaxBinVAll = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
333  m_stripMaxBinVAll->GetXaxis()->SetTitle("max bin");
334  m_stripMaxBinVAll->GetYaxis()->SetTitle("count");
335  m_histoList->Add(m_stripMaxBinVAll);
336 
337  name = str(format("SVDDQM_StripMaxBinU3"));
338  title = str(format("SVD U-Strip MaxBin for layer 3 sensors"));
339  m_stripMaxBinU3 = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
340  m_stripMaxBinU3->GetXaxis()->SetTitle("max bin");
341  m_stripMaxBinU3->GetYaxis()->SetTitle("count");
342  m_histoList->Add(m_stripMaxBinU3);
343  name = str(format("SVDDQM_StripMaxBinV3"));
344  title = str(format("SVD V-Strip MaxBin for layer 3 sensors"));
345  m_stripMaxBinV3 = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
346  m_stripMaxBinV3->GetXaxis()->SetTitle("max bin");
347  m_stripMaxBinV3->GetYaxis()->SetTitle("count");
348  m_histoList->Add(m_stripMaxBinV3);
349 
350  name = str(format("SVDDQM_StripMaxBinU6"));
351  title = str(format("SVD U-Strip MaxBin for layer 6 sensors"));
352  m_stripMaxBinU6 = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
353  m_stripMaxBinU6->GetXaxis()->SetTitle("max bin");
354  m_stripMaxBinU6->GetYaxis()->SetTitle("count");
355  m_histoList->Add(m_stripMaxBinU6);
356  name = str(format("SVDDQM_StripMaxBinV6"));
357  title = str(format("SVD V-Strip MaxBin for layer 6 sensors"));
358  m_stripMaxBinV6 = new TH1F(name.c_str(), title.c_str(), MaxBinBins, 0, MaxBinMax);
359  m_stripMaxBinV6->GetXaxis()->SetTitle("max bin");
360  m_stripMaxBinV6->GetYaxis()->SetTitle("count");
361  m_histoList->Add(m_stripMaxBinV6);
362 
363  for (int i = 0; i < nSVDSensors; i++) {
364  VxdID id = gTools->getSensorIDFromSVDIndex(i);
365  int iLayer = id.getLayerNumber();
366  int iLadder = id.getLadderNumber();
367  int iSensor = id.getSensorNumber();
368  VxdID sensorID(iLayer, iLadder, iSensor);
369  SVD::SensorInfo SensorInfo = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::get(sensorID));
370  string sensorDescr = str(format("%1%_%2%_%3%") % iLayer % iLadder % iSensor);
371  //----------------------------------------------------------------
372  // Number of fired strips per sensor
373  //----------------------------------------------------------------
374  name = str(format("SVDDQM_%1%_FiredU") % sensorDescr);
375  title = str(format("SVD Sensor %1% Number of Fired U-Strips") % sensorDescr);
376  m_firedU[i] = new TH1F(name.c_str(), title.c_str(), 50, 0, 50);
377  m_firedU[i]->GetXaxis()->SetTitle("# fired strips");
378  m_firedU[i]->GetYaxis()->SetTitle("count");
379  m_histoList->Add(m_firedU[i]);
380  name = str(format("SVDDQM_%1%_FiredV") % sensorDescr);
381  title = str(format("SVD Sensor %1% Number of Fired V-Strips") % sensorDescr);
382  m_firedV[i] = new TH1F(name.c_str(), title.c_str(), 50, 0, 50);
383  m_firedV[i]->GetXaxis()->SetTitle("# fired strips");
384  m_firedV[i]->GetYaxis()->SetTitle("count");
385  m_histoList->Add(m_firedV[i]);
386  //----------------------------------------------------------------
387  // Number of clusters per sensor
388  //----------------------------------------------------------------
389  name = str(format("SVDDQM_%1%_ClustersU") % sensorDescr);
390  title = str(format("SVD Sensor %1% Number of U-Clusters") % sensorDescr);
391  m_clustersU[i] = new TH1F(name.c_str(), title.c_str(), 20, 0, 20);
392  m_clustersU[i]->GetXaxis()->SetTitle("# clusters");
393  m_clustersU[i]->GetYaxis()->SetTitle("count");
394  m_histoList->Add(m_clustersU[i]);
395  name = str(format("SVDDQM_%1%_ClustersV") % sensorDescr);
396  title = str(format("SVD Sensor %1% Number of V-Clusters") % sensorDescr);
397  m_clustersV[i] = new TH1F(name.c_str(), title.c_str(), 20, 0, 20);
398  m_clustersV[i]->GetXaxis()->SetTitle("# clusters");
399  m_clustersV[i]->GetYaxis()->SetTitle("count");
400  m_histoList->Add(m_clustersV[i]);
401  //----------------------------------------------------------------
402  // Charge of clusters
403  //----------------------------------------------------------------
404  name = str(format("SVDDQM_%1%_ClusterChargeU") % sensorDescr);
405  title = str(format("SVD Sensor %1% U-Cluster Charge") % sensorDescr);
406  m_clusterChargeU[i] = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
407  m_clusterChargeU[i]->GetXaxis()->SetTitle("cluster charge [ke-]");
408  m_clusterChargeU[i]->GetYaxis()->SetTitle("count");
409  m_histoList->Add(m_clusterChargeU[i]);
410  name = str(format("SVDDQM_%1%_ClusterChargeV") % sensorDescr);
411  title = str(format("SVD Sensor %1% V-Cluster Charge") % sensorDescr);
412  m_clusterChargeV[i] = new TH1F(name.c_str(), title.c_str(), ChargeBins, 0, ChargeMax);
413  m_clusterChargeV[i]->GetXaxis()->SetTitle("cluster charge [ke-]");
414  m_clusterChargeV[i]->GetYaxis()->SetTitle("count");
415  m_histoList->Add(m_clusterChargeV[i]);
416  //----------------------------------------------------------------
417  // SNR of clusters
418  //----------------------------------------------------------------
419  name = str(format("SVDDQM_%1%_ClusterSNRU") % sensorDescr);
420  title = str(format("SVD Sensor %1% U-Cluster SNR") % sensorDescr);
421  m_clusterSNRU[i] = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
422  m_clusterSNRU[i]->GetXaxis()->SetTitle("cluster SNR");
423  m_clusterSNRU[i]->GetYaxis()->SetTitle("count");
424  m_histoList->Add(m_clusterSNRU[i]);
425  name = str(format("SVDDQM_%1%_ClusterSNRV") % sensorDescr);
426  title = str(format("SVD Sensor %1% V-Cluster SNR") % sensorDescr);
427  m_clusterSNRV[i] = new TH1F(name.c_str(), title.c_str(), SNRBins, 0, SNRMax);
428  m_clusterSNRV[i]->GetXaxis()->SetTitle("cluster SNR");
429  m_clusterSNRV[i]->GetYaxis()->SetTitle("count");
430  m_histoList->Add(m_clusterSNRV[i]);
431  //----------------------------------------------------------------
432  // Charge of strips
433  //----------------------------------------------------------------
434  name = str(format("SVDDQM_%1%_ADCStripU") % sensorDescr);
435  title = str(format("SVD Sensor %1% U-Strip signal in ADC Counts, all 6 APV samples") % sensorDescr);
436  m_stripSignalU[i] = new TH1F(name.c_str(), title.c_str(), 256, -0.5, 255.5);
437  m_stripSignalU[i]->GetXaxis()->SetTitle("signal ADC");
438  m_stripSignalU[i]->GetYaxis()->SetTitle("count");
439  m_histoList->Add(m_stripSignalU[i]);
440  name = str(format("SVDDQM_%1%_ADCStripV") % sensorDescr);
441  title = str(format("SVD Sensor %1% V-Strip signal in ADC Counts, all 6 APV samples") % sensorDescr);
442  m_stripSignalV[i] = new TH1F(name.c_str(), title.c_str(), 256, -0.5, 255.5);
443  m_stripSignalV[i]->GetXaxis()->SetTitle("signal ADC");
444  m_stripSignalV[i]->GetYaxis()->SetTitle("count");
445  m_histoList->Add(m_stripSignalV[i]);
446  //----------------------------------------------------------------
447  // Strips Counts
448  //----------------------------------------------------------------
449  name = str(format("SVDDQM_%1%_StripCountU") % sensorDescr);
450  title = str(format("SVD Sensor %1% Integrated Number of ZS5 Fired U-Strip vs Strip Number") % sensorDescr);
451  m_stripCountU[i] = new TH1F(name.c_str(), title.c_str(), 768, -0.5, 767.5);
452  m_stripCountU[i]->GetXaxis()->SetTitle("cellID");
453  m_stripCountU[i]->GetYaxis()->SetTitle("count");
454  m_histoList->Add(m_stripCountU[i]);
455  name = str(format("SVDDQM_%1%_StripCountV") % sensorDescr);
456  title = str(format("SVD Sensor %1% Integrated Number of ZS5 Fired V-Strip vs Strip Number") % sensorDescr);
457  m_stripCountV[i] = new TH1F(name.c_str(), title.c_str(), 768, -0.5, 767.5);
458  m_stripCountV[i]->GetXaxis()->SetTitle("cellID");
459  m_stripCountV[i]->GetYaxis()->SetTitle("count");
460  m_histoList->Add(m_stripCountV[i]);
461  //----------------------------------------------------------------
462  // Strips Counts with online ZS
463  //----------------------------------------------------------------
464  name = str(format("SVDDQM_%1%_OnlineZSStripCountU") % sensorDescr);
465  title = str(format("SVD Sensor %1% Integrated Number of online-ZS Fired U-Strip vs Strip Number") % sensorDescr);
466  m_onlineZSstripCountU[i] = new TH1F(name.c_str(), title.c_str(), 768, -0.5, 767.5);
467  m_onlineZSstripCountU[i]->GetXaxis()->SetTitle("cellID");
468  m_onlineZSstripCountU[i]->GetYaxis()->SetTitle("count");
469  m_histoList->Add(m_onlineZSstripCountU[i]);
470  name = str(format("SVDDQM_%1%_OnlineZSStripCountV") % sensorDescr);
471  title = str(format("SVD Sensor %1% Integrated Number of online-ZS Fired V-Strip vs Strip Number") % sensorDescr);
472  m_onlineZSstripCountV[i] = new TH1F(name.c_str(), title.c_str(), 768, -0.5, 767.5);
473  m_onlineZSstripCountV[i]->GetXaxis()->SetTitle("cellID");
474  m_onlineZSstripCountV[i]->GetYaxis()->SetTitle("count");
475  m_histoList->Add(m_onlineZSstripCountV[i]);
476  //----------------------------------------------------------------
477  // Cluster size distribution
478  //----------------------------------------------------------------
479  name = str(format("SVDDQM_%1%_ClusterSizeU") % sensorDescr);
480  title = str(format("SVD Sensor %1% U-Cluster Size") % sensorDescr);
481  m_clusterSizeU[i] = new TH1F(name.c_str(), title.c_str(), 9, 1, 10);
482  m_clusterSizeU[i]->GetXaxis()->SetTitle("cluster size");
483  m_clusterSizeU[i]->GetYaxis()->SetTitle("count");
484  m_histoList->Add(m_clusterSizeU[i]);
485  name = str(format("SVDDQM_%1%_ClusterSizeV") % sensorDescr);
486  title = str(format("SVD Sensor %1% V-Cluster Size") % sensorDescr);
487  m_clusterSizeV[i] = new TH1F(name.c_str(), title.c_str(), 9, 1, 10);
488  m_clusterSizeV[i]->GetXaxis()->SetTitle("cluster size");
489  m_clusterSizeV[i]->GetYaxis()->SetTitle("count");
490  m_histoList->Add(m_clusterSizeV[i]);
491  //----------------------------------------------------------------
492  // Cluster time distribution
493  //----------------------------------------------------------------
494  name = str(format("SVDDQM_%1%_ClusterTimeU") % sensorDescr);
495  Title = Form("SVD Sensor %s U-Cluster Time %s", sensorDescr.c_str(), refFrame.Data());
496  m_clusterTimeU[i] = new TH1F(name.c_str(), Title.Data(), TimeBins, TimeMin, TimeMax);
497  m_clusterTimeU[i]->GetXaxis()->SetTitle("cluster time (ns)");
498  m_clusterTimeU[i]->GetYaxis()->SetTitle("count");
499  m_histoList->Add(m_clusterTimeU[i]);
500  name = str(format("SVDDQM_%1%_ClusterTimeV") % sensorDescr);
501  Title = Form("SVD Sensor %s V-Cluster Time %s", sensorDescr.c_str(), refFrame.Data());
502  m_clusterTimeV[i] = new TH1F(name.c_str(), Title.Data(), TimeBins, TimeMin, TimeMax);
503  m_clusterTimeV[i]->GetXaxis()->SetTitle("cluster time (ns)");
504  m_clusterTimeV[i]->GetYaxis()->SetTitle("count");
505  m_histoList->Add(m_clusterTimeV[i]);
506  }
507 
508  for (int i = 0; i < nSVDChips; i++) {
509  VxdID id = gTools->getChipIDFromSVDIndex(i);
510  int iLayer = id.getLayerNumber();
511  int iLadder = id.getLadderNumber();
512  int iSensor = id.getSensorNumber();
513  int iChip = gTools->getSVDChipNumber(id);
514  int IsU = gTools->isSVDSideU(id);
515  TString AxisTicks = Form("%i_%i_%i_u%i", iLayer, iLadder, iSensor, iChip);
516  if (!IsU)
517  AxisTicks = Form("%i_%i_%i_v%i", iLayer, iLadder, iSensor, iChip);
518  m_hitMapCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
519  m_hitMapClCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
520  }
521 
522 
523 
524  //----------------------------------------------------------------
525  // Additional histograms for out of ExpressReco
526  //----------------------------------------------------------------
527 
528  if (m_ShowAllHistos == 1) {
529  TDirectory* dirShowAll = NULL;
530  dirShowAll = oldDir->mkdir("SVDDQMAll");
531  dirShowAll->cd();
532 
533  m_hitMapU = new TH2F*[nSVDSensors];
534  m_hitMapV = new TH2F*[nSVDSensors];
535  m_hitMapUCl = new TH1F*[nSVDSensors];
536  m_hitMapVCl = new TH1F*[nSVDSensors];
537  for (int i = 0; i < nSVDSensors; i++) {
538  VxdID id = gTools->getSensorIDFromSVDIndex(i);
539  int iLayer = id.getLayerNumber();
540  int iLadder = id.getLadderNumber();
541  int iSensor = id.getSensorNumber();
542  VxdID sensorID(iLayer, iLadder, iSensor);
543  SVD::SensorInfo SensorInfo = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::get(sensorID));
544  string sensorDescr = str(format("%1%_%2%_%3%") % iLayer % iLadder % iSensor);
545  //----------------------------------------------------------------
546  // Hitmaps: Number of strips by coordinate
547  //----------------------------------------------------------------
548  name = str(format("SVD_%1%_StripHitmapU") % sensorDescr);
549  title = str(format("SVD Sensor %1% Strip Hitmap in U") % sensorDescr);
550  int nStrips = SensorInfo.getUCells();
551  m_hitMapU[i] = new TH2F(name.c_str(), title.c_str(), nStrips, 0, nStrips, SVDShaperDigit::c_nAPVSamples, 0,
552  SVDShaperDigit::c_nAPVSamples);
553  m_hitMapU[i]->GetXaxis()->SetTitle("u position [pitch units]");
554  m_hitMapU[i]->GetYaxis()->SetTitle("timebin [time units]");
555  m_hitMapU[i]->GetZaxis()->SetTitle("hits");
556  m_histoList->Add(m_hitMapU[i]);
557  name = str(format("SVD_%1%_StripHitmapV") % sensorDescr);
558  title = str(format("SVD Sensor %1% Strip Hitmap in V") % sensorDescr);
559  nStrips = SensorInfo.getVCells();
560  m_hitMapV[i] = new TH2F(name.c_str(), title.c_str(), nStrips, 0, nStrips, SVDShaperDigit::c_nAPVSamples, 0,
561  SVDShaperDigit::c_nAPVSamples);
562  m_hitMapV[i]->GetXaxis()->SetTitle("v position [pitch units]");
563  m_hitMapV[i]->GetYaxis()->SetTitle("timebin [time units]");
564  m_hitMapV[i]->GetZaxis()->SetTitle("hits");
565  m_histoList->Add(m_hitMapV[i]);
566  //----------------------------------------------------------------
567  // Hitmaps: Number of clusters by coordinate
568  //----------------------------------------------------------------
569  name = str(format("SVD_%1%_HitmapClstU") % sensorDescr);
570  title = str(format("SVD Sensor %1% Hitmap Clusters in U") % sensorDescr);
571  nStrips = SensorInfo.getUCells();
572  m_hitMapUCl[i] = new TH1F(name.c_str(), title.c_str(), nStrips, 0, nStrips);
573  m_hitMapUCl[i]->GetXaxis()->SetTitle("u position [pitch units]");
574  m_hitMapUCl[i]->GetYaxis()->SetTitle("hits");
575  m_histoList->Add(m_hitMapUCl[i]);
576  name = str(format("SVD_%1%_HitmapClstV") % sensorDescr);
577  title = str(format("SVD Sensor %1% Hitmap Clusters in V") % sensorDescr);
578  nStrips = SensorInfo.getVCells();
579  m_hitMapVCl[i] = new TH1F(name.c_str(), title.c_str(), nStrips, 0, nStrips);
580  m_hitMapVCl[i]->GetXaxis()->SetTitle("v position [pitch units]");
581  m_hitMapVCl[i]->GetYaxis()->SetTitle("hits");
582  m_histoList->Add(m_hitMapVCl[i]);
583  }
584  }
585 
586  oldDir->cd();
587 }
588 
589 
590 void SVDDQMExpressRecoModule::initialize()
591 {
592  // Register histograms (calls back defineHisto)
593  REG_HISTOGRAM
594 
595  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
596  if (gTools->getNumberOfSVDLayers() != 0) {
597  //Register collections
598  StoreArray<SVDShaperDigit> storeNoZSSVDShaperDigits(m_storeNoZSSVDShaperDigitsName);
599  StoreArray<SVDShaperDigit> storeSVDShaperDigits(m_storeSVDShaperDigitsName);
600  StoreArray<SVDCluster> storeSVDClusters(m_storeSVDClustersName);
601  m_storeSVDClustersName = storeSVDClusters.getName();
602 
603  storeSVDClusters.isOptional();
604  storeSVDShaperDigits.isOptional();
605  m_svdEventInfo.isOptional();
606  storeNoZSSVDShaperDigits.isOptional();
607 
608  //Store names to speed up creation later
609  m_storeSVDShaperDigitsName = storeSVDShaperDigits.getName();
610  }
611 }
612 
613 void SVDDQMExpressRecoModule::beginRun()
614 {
615  StoreObjPtr<EventMetaData> evtMetaData;
616  m_expNumber = evtMetaData->getExperiment();
617  m_runNumber = evtMetaData->getRun();
618 
619  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
620  if (gTools->getNumberOfSVDLayers() == 0) return;
621 
622  // Add experiment and run number to the title of selected histograms (CR shifter plots)
623  TString runID = TString::Format(" ~ Exp%d Run%d", m_expNumber, m_runNumber);
624  TObject* obj;
625  TIter nextH(m_histoList);
626  while ((obj = nextH()))
627  if (obj->InheritsFrom("TH1")) {
628  ((TH1F*)obj)->SetTitle(obj->GetTitle() + runID);
629  if (obj != NULL)((TH1F*)obj)->Reset();
630  }
631 
632 }
633 
634 void SVDDQMExpressRecoModule::event()
635 {
636 
637  //check HLT decision and increase number of events only if the event has been accepted
638 
639  if (m_skipRejectedEvents && (m_resultStoreObjectPointer.isValid())) {
640  const bool eventAccepted = FinalTriggerDecisionCalculator::getFinalTriggerDecision(*m_resultStoreObjectPointer);
641  if (!eventAccepted) return;
642  }
643  m_nEvents->Fill(0);
644 
645  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
646  if (gTools->getNumberOfSVDLayers() == 0) return;
647 
648  const StoreArray<SVDShaperDigit> storeNoZSSVDShaperDigits(m_storeNoZSSVDShaperDigitsName);
649 
650  const StoreArray<SVDShaperDigit> storeSVDShaperDigits(m_storeSVDShaperDigitsName);
651  const StoreArray<SVDCluster> storeSVDClusters(m_storeSVDClustersName);
652 
653  if (!storeSVDShaperDigits.isValid() || !storeSVDShaperDigits.getEntries()) {
654  return;
655  }
656 
657  int firstSVDLayer = gTools->getFirstSVDLayer();
658  int lastSVDLayer = gTools->getLastSVDLayer();
659  int nSVDSensors = gTools->getNumberOfSVDSensors();
660 
661  // Fired strips offline ZS
662  vector< set<int> > uStrips(nSVDSensors); // sets to eliminate multiple samples per strip
663  vector< set<int> > vStrips(nSVDSensors);
664  for (const SVDShaperDigit& digitIn : storeSVDShaperDigits) {
665  int iLayer = digitIn.getSensorID().getLayerNumber();
666  if ((iLayer < firstSVDLayer) || (iLayer > lastSVDLayer)) continue;
667  int iLadder = digitIn.getSensorID().getLadderNumber();
668  int iSensor = digitIn.getSensorID().getSensorNumber();
669  VxdID sensorID(iLayer, iLadder, iSensor);
670  int index = gTools->getSVDSensorIndex(sensorID);
671  SVD::SensorInfo SensorInfo = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::get(sensorID));
672  if (digitIn.isUStrip()) {
673 
674  // B2DEBUG(29, digitIn.toString().c_str() );
675 
676  //fill strip count first
677  if (m_stripCountU[index] != NULL) m_stripCountU[index]->Fill(digitIn.getCellID());
678 
679  //fill max bin
680  if (m_stripMaxBinUAll != NULL) m_stripMaxBinUAll->Fill(digitIn.getMaxTimeBin());
681  if (iLayer == 3)
682  if (m_stripMaxBinU3 != NULL) m_stripMaxBinU3->Fill(digitIn.getMaxTimeBin());
683  if (iLayer == 6)
684  if (m_stripMaxBinU6 != NULL) m_stripMaxBinU6->Fill(digitIn.getMaxTimeBin());
685 
686  uStrips.at(index).insert(digitIn.getCellID());
687  int Chip = (int)(digitIn.getCellID() / gTools->getSVDChannelsPerChip()) + 1;
688  int indexChip = gTools->getSVDChipIndex(sensorID, kTRUE, Chip);
689  // 6-to-1 relation weights are equal to digit signals, modulo rounding error
690  SVDShaperDigit::APVFloatSamples samples = digitIn.getSamples();
691  int isSample = 0;
692  for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
693  if (m_stripSignalU[index] != NULL) m_stripSignalU[index]->Fill(samples[i]);
694  if (samples[i] > m_CutSVDCharge) {
695  isSample = 1;
696  if (m_ShowAllHistos == 1) {
697  if (m_hitMapU[index] != NULL) m_hitMapU[index]->Fill(digitIn.getCellID(), i);
698  }
699  }
700  }
701  if (isSample) {
702  if (m_hitMapCountsU != NULL) m_hitMapCountsU->Fill(index);
703  if (m_hitMapCountsChip != NULL) m_hitMapCountsChip->Fill(indexChip);
704  }
705  } else {
706  //fill strip count first
707  if (m_stripCountV[index] != NULL) m_stripCountV[index]->Fill(digitIn.getCellID());
708 
709  //fill max bin
710  if (m_stripMaxBinVAll != NULL) m_stripMaxBinVAll->Fill(digitIn.getMaxTimeBin());
711 
712  if (iLayer == 3)
713  if (m_stripMaxBinV3 != NULL) m_stripMaxBinV3->Fill(digitIn.getMaxTimeBin());
714  if (iLayer == 6)
715  if (m_stripMaxBinV6 != NULL) m_stripMaxBinV6->Fill(digitIn.getMaxTimeBin());
716 
717  vStrips.at(index).insert(digitIn.getCellID());
718  int Chip = (int)(digitIn.getCellID() / gTools->getSVDChannelsPerChip()) + 1;
719  int indexChip = gTools->getSVDChipIndex(sensorID, kFALSE, Chip);
720  // 6-to-1 relation weights are equal to digit signals, modulo rounding error
721  SVDShaperDigit::APVFloatSamples samples = digitIn.getSamples();
722  int isSample = 0;
723  for (size_t i = 0; i < SVDShaperDigit::c_nAPVSamples; ++i) {
724  if (m_stripSignalV[index] != NULL) m_stripSignalV[index]->Fill(samples[i]);
725  if (samples[i] > m_CutSVDCharge) {
726  isSample = 1;
727  if (m_ShowAllHistos == 1) {
728  if (m_hitMapV[index] != NULL) m_hitMapV[index]->Fill(digitIn.getCellID(), i);
729  }
730  }
731  }
732  if (isSample) {
733  if (m_hitMapCountsV != NULL) m_hitMapCountsV->Fill(index);
734  if (m_hitMapCountsChip != NULL) m_hitMapCountsChip->Fill(indexChip);
735  }
736  }
737  }
738  for (int i = 0; i < nSVDSensors; i++) {
739  if ((m_firedU[i] != NULL) && (uStrips[i].size() > 0))
740  m_firedU[i]->Fill(uStrips[i].size());
741  if ((m_firedV[i] != NULL) && (vStrips[i].size() > 0))
742  m_firedV[i]->Fill(vStrips[i].size());
743  }
744 
745  // Fired strips ONLINE ZS
746  if (storeNoZSSVDShaperDigits.isValid())
747  for (const SVDShaperDigit& digitIn : storeNoZSSVDShaperDigits) {
748  int iLayer = digitIn.getSensorID().getLayerNumber();
749  if ((iLayer < firstSVDLayer) || (iLayer > lastSVDLayer)) continue;
750  int iLadder = digitIn.getSensorID().getLadderNumber();
751  int iSensor = digitIn.getSensorID().getSensorNumber();
752  VxdID sensorID(iLayer, iLadder, iSensor);
753  int index = gTools->getSVDSensorIndex(sensorID);
754  SVD::SensorInfo SensorInfo = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::get(sensorID));
755  if (digitIn.isUStrip()) {
756  if (m_onlineZSstripCountU[index] != NULL) m_onlineZSstripCountU[index]->Fill(digitIn.getCellID());
757  } else {
758  if (m_onlineZSstripCountV[index] != NULL) m_onlineZSstripCountV[index]->Fill(digitIn.getCellID());
759  }
760  }
761 
762  vector< set<int> > countsU(nSVDSensors); // sets to eliminate multiple samples per strip
763  vector< set<int> > countsV(nSVDSensors);
764  // Hitmaps, Charge, Seed, Size, Time, ...
765  for (const SVDCluster& cluster : storeSVDClusters) {
766  if (cluster.getCharge() < m_CutSVDClusterCharge) continue;
767  int iLayer = cluster.getSensorID().getLayerNumber();
768  if ((iLayer < firstSVDLayer) || (iLayer > lastSVDLayer)) continue;
769  int iLadder = cluster.getSensorID().getLadderNumber();
770  int iSensor = cluster.getSensorID().getSensorNumber();
771  VxdID sensorID(iLayer, iLadder, iSensor);
772  int index = gTools->getSVDSensorIndex(sensorID);
773  SVD::SensorInfo SensorInfo = dynamic_cast<const SVD::SensorInfo&>(VXD::GeoCache::get(sensorID));
774 
775  float time = cluster.getClsTime();
776  if (m_desynchSVDTime && m_svdEventInfo.isValid())
777  time = time - m_svdEventInfo->getSVD2FTSWTimeShift(cluster.getFirstFrame());
778 
779  if (cluster.isUCluster()) {
780  countsU.at(index).insert(SensorInfo.getUCellID(cluster.getPosition()));
781  int indexChip = gTools->getSVDChipIndex(sensorID, kTRUE,
782  (int)(SensorInfo.getUCellID(cluster.getPosition()) / gTools->getSVDChannelsPerChip()) + 1);
783  if (m_hitMapClCountsU != NULL) m_hitMapClCountsU->Fill(index);
784  if (m_hitMapClCountsChip != NULL) m_hitMapClCountsChip->Fill(indexChip);
785  if (m_clusterChargeU[index] != NULL) m_clusterChargeU[index]->Fill(cluster.getCharge() / 1000.0); // in kelectrons
786  if (m_clusterSNRU[index] != NULL) m_clusterSNRU[index]->Fill(cluster.getSNR());
787  if (m_clusterChargeUAll != NULL) m_clusterChargeUAll->Fill(cluster.getCharge() / 1000.0); // in kelectrons
788  if (m_clusterSNRUAll != NULL) m_clusterSNRUAll->Fill(cluster.getSNR());
789  if (m_clusterSizeU[index] != NULL) m_clusterSizeU[index]->Fill(cluster.getSize());
790  if (m_clusterTimeU[index] != NULL) m_clusterTimeU[index]->Fill(time);
791  if (m_clusterTimeUAll != NULL) m_clusterTimeUAll->Fill(time);
792  if (iLayer == 3) {
793  if (m_clusterChargeU3 != NULL) m_clusterChargeU3->Fill(cluster.getCharge() / 1000.0); // in kelectrons
794  if (m_clusterSNRU3 != NULL) m_clusterSNRU3->Fill(cluster.getSNR());
795  if (m_clusterTimeU3 != NULL) m_clusterTimeU3->Fill(time);
796  } else {
797  if (m_clusterChargeU456 != NULL) m_clusterChargeU456->Fill(cluster.getCharge() / 1000.0); // in kelectrons
798  if (m_clusterSNRU456 != NULL) m_clusterSNRU456->Fill(cluster.getSNR());
799  if (m_clusterTimeU456 != NULL) m_clusterTimeU456->Fill(time);
800  }
801 
802  if (m_ShowAllHistos == 1)
803  if (m_hitMapUCl[index] != NULL) m_hitMapUCl[index]->Fill(SensorInfo.getUCellID(cluster.getPosition()));
804  } else {
805  countsV.at(index).insert(SensorInfo.getVCellID(cluster.getPosition()));
806  int indexChip = gTools->getSVDChipIndex(sensorID, kFALSE,
807  (int)(SensorInfo.getVCellID(cluster.getPosition()) / gTools->getSVDChannelsPerChip()) + 1);
808  if (m_hitMapClCountsV != NULL) m_hitMapClCountsV->Fill(index);
809  if (m_hitMapClCountsChip != NULL) m_hitMapClCountsChip->Fill(indexChip);
810  if (m_clusterChargeV[index] != NULL) m_clusterChargeV[index]->Fill(cluster.getCharge() / 1000.0); // in kelectrons
811  if (m_clusterSNRV[index] != NULL) m_clusterSNRV[index]->Fill(cluster.getSNR());
812  if (m_clusterChargeVAll != NULL) m_clusterChargeVAll->Fill(cluster.getCharge() / 1000.0); // in kelectrons
813  if (m_clusterSNRVAll != NULL) m_clusterSNRVAll->Fill(cluster.getSNR());
814  if (m_clusterSizeV[index] != NULL) m_clusterSizeV[index]->Fill(cluster.getSize());
815  if (m_clusterTimeV[index] != NULL) m_clusterTimeV[index]->Fill(time);
816  if (m_clusterTimeVAll != NULL) m_clusterTimeVAll->Fill(time);
817  if (iLayer == 3) {
818  if (m_clusterChargeV3 != NULL) m_clusterChargeV3->Fill(cluster.getCharge() / 1000.0); // in kelectrons
819  if (m_clusterSNRV3 != NULL) m_clusterSNRV3->Fill(cluster.getSNR());
820  if (m_clusterTimeV3 != NULL) m_clusterTimeV3->Fill(time);
821  } else {
822  if (m_clusterChargeV456 != NULL) m_clusterChargeV456->Fill(cluster.getCharge() / 1000.0); // in kelectrons
823  if (m_clusterSNRV456 != NULL) m_clusterSNRV456->Fill(cluster.getSNR());
824  if (m_clusterTimeV456 != NULL) m_clusterTimeV456->Fill(time);
825  }
826  if (m_ShowAllHistos == 1)
827  if (m_hitMapVCl[index] != NULL) m_hitMapVCl[index]->Fill(SensorInfo.getVCellID(cluster.getPosition()));
828 
829  }
830  }
831  for (int i = 0; i < nSVDSensors; i++) {
832  if ((m_clustersU[i] != NULL) && (countsU[i].size() > 0))
833  m_clustersU[i]->Fill(countsU[i].size());
834  if ((m_clustersV[i] != NULL) && (countsV[i].size() > 0))
835  m_clustersV[i]->Fill(countsV[i].size());
836  }
837 }
838 
839 
840 void SVDDQMExpressRecoModule::terminate()
841 {
842 
843  // m_histoList->Delete();
844  delete m_histoList;
845 
846 }
Belle2::VXD::SensorInfoBase::getUCells
int getUCells() const
Return number of pixel/strips in u direction.
Definition: SensorInfoBase.h:223
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::VXD::SensorInfoBase::getUCellID
int getUCellID(double u, double v=0, bool clamp=false) const
Return the corresponding pixel/strip ID of a given u coordinate.
Definition: SensorInfoBase.h:203
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::SVD::SensorInfo
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
Definition: SensorInfo.h:35
Belle2::SVDShaperDigit
The SVD ShaperDigit class.
Definition: SVDShaperDigit.h:46
Belle2::SVDShaperDigit::APVFloatSamples
std::array< APVFloatSampleType, c_nAPVSamples > APVFloatSamples
array of APVFloatSampleType objects
Definition: SVDShaperDigit.h:63
Belle2::PXD::SensorInfo
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:34
Belle2::VXD::SensorInfoBase::getVCellID
int getVCellID(double v, bool clamp=false) const
Return the corresponding pixel/strip ID of a given v coordinate.
Definition: SensorInfoBase.h:213
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::VXD::SensorInfoBase::getVCells
int getVCells() const
Return number of pixel/strips in v direction.
Definition: SensorInfoBase.h:225
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::StoreArray::isValid
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:298
Belle2::SVDDQMExpressRecoModule
SVD DQM Module for Express Reco.
Definition: SVDDQMExpressRecoModule.h:47
Belle2::SVDCluster
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:38
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29