Belle II Software  release-05-02-19
PXDDQMClustersModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Björn Spruck, *
7  * *
8  * Prepared for Belle II geometry *
9  * *
10  * This software is provided "as is" without any warranty. *
11  **************************************************************************/
12 
13 #include "pxd/modules/pxdDQM/PXDDQMClustersModule.h"
14 
15 #include <framework/core/HistoModule.h>
16 #include <framework/gearbox/Unit.h>
17 #include <framework/datastore/StoreObjPtr.h>
18 #include <framework/datastore/StoreArray.h>
19 #include <framework/datastore/RelationArray.h>
20 
21 #include <pxd/dataobjects/PXDDigit.h>
22 #include <pxd/dataobjects/PXDCluster.h>
23 #include <pxd/geometry/SensorInfo.h>
24 
25 #include <vxd/geometry/GeoCache.h>
26 #include <vxd/geometry/SensorInfoBase.h>
27 #include <vxd/geometry/GeoTools.h>
28 #include <pxd/unpacking/PXDMappingLookup.h>
29 #include <pxd/reconstruction/PXDGainCalibrator.h>
30 
31 #include <boost/format.hpp>
32 
33 #include "TDirectory.h"
34 
35 using namespace std;
36 using boost::format;
37 using namespace Belle2;
38 using namespace Belle2::PXD;
39 
40 //-----------------------------------------------------------------
41 // Register the Module
42 //-----------------------------------------------------------------
43 REG_MODULE(PXDDQMClusters)
44 
45 
46 //-----------------------------------------------------------------
47 // Implementation
48 //-----------------------------------------------------------------
49 
51 {
52  //Set module properties
53  setDescription("PXD DQM clusters module "
54  "Recommended Number of events for monitorin is 40 kEvents or more to fill all histograms "
55  );
56 
57  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
58  addParam("CutMinCharge", m_CutMinCharge,
59  "cut on pixel charge for accepting good pixel, default >= 12", 12);
60  addParam("CutMinSeedCharge", m_CutMinSeedCharge,
61  "cut on cluster seed for accepting good cluster, default >= 12", 12);
62  addParam("CutMaxClusterSize", m_CutMaxClusterSize,
63  "cut on cluster size accepting good cluster, default <= 4", 4);
64  addParam("CutMinClusterCharge", m_CutMinClusterCharge,
65  "cut on cluster charge for accepting good cluster, default >= 12", 12);
66  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
67  std::string("PXDDQMClusters"));
68 }
69 
70 
71 //------------------------------------------------------------------
72 // Function to define histograms
73 //-----------------------------------------------------------------
74 
75 void PXDDQMClustersModule::defineHisto()
76 {
77  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
78  if (gTools->getNumberOfLayers() == 0) {
79  B2FATAL("Missing geometry for VXD, check steering file.");
80  }
81  if (gTools->getNumberOfPXDLayers() == 0) {
82  B2WARNING("Missing geometry for PXD, PXD-DQM is skipped.");
83  return;
84  }
85 
86  // Create a separate histogram directories and cd into it.
87  TDirectory* oldDir = gDirectory;
88  if (m_histogramDirectoryName != "") {
89  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
90  oldDir->cd(m_histogramDirectoryName.c_str());
91  }
92 
93  // basic constants presets:
94  int nPXDSensors = gTools->getNumberOfPXDSensors();
95  int nPXDChips = gTools->getTotalPXDChips();
96 
97  // Create basic histograms:
98  m_hitMapCounts = new TH1D("DQM_PXD_PixelHitmapCounts", "PXD Integrated number of fired pixels per sensor",
99  nPXDSensors, 0, nPXDSensors);
100  m_hitMapCounts->GetXaxis()->SetTitle("Sensor ID");
101  m_hitMapCounts->GetYaxis()->SetTitle("counts");
102 
103  m_hitMapFilterCounts = new TH1D("DQM_PXD_PixelHitmapFilterCounts", "PXD Integrated number of filtered pixels per sensor",
104  nPXDSensors, 0, nPXDSensors);
105  m_hitMapFilterCounts->GetXaxis()->SetTitle("Sensor ID");
106  m_hitMapFilterCounts->GetYaxis()->SetTitle("counts");
107 
108  m_hitMapClCounts = new TH1D("DQM_PXD_ClusterHitmapCounts", "PXD Integrated number of clusters per sensor",
109  nPXDSensors, 0, nPXDSensors);
110  m_hitMapClCounts->GetXaxis()->SetTitle("Sensor ID");
111  m_hitMapClCounts->GetYaxis()->SetTitle("counts");
112 
113  m_hitMapClFilterCounts = new TH1D("DQM_PXD_ClusterHitmapFilterCounts", "PXD Integrated number of filtered clusters per sensor",
114  nPXDSensors, 0, nPXDSensors);
115  m_hitMapClFilterCounts->GetXaxis()->SetTitle("Sensor ID");
116  m_hitMapClFilterCounts->GetYaxis()->SetTitle("counts");
117 
118  // basic counters per chip:
119  m_hitMapCountsChip = new TH1D("DQM_PXD_PixelHitmapCountsChip", "PXD Integrated number of fired pixels per chip",
120  nPXDChips, 0, nPXDChips);
121  m_hitMapCountsChip->GetXaxis()->SetTitle("Chip ID");
122  m_hitMapCountsChip->GetYaxis()->SetTitle("counts");
123  m_hitMapClCountsChip = new TH1D("DQM_PXD_ClusterHitmapCountsChip", "PXD Integrated number of clusters per chip",
124  nPXDChips, 0, nPXDChips);
125  m_hitMapClCountsChip->GetXaxis()->SetTitle("Chip ID");
126  m_hitMapClCountsChip->GetYaxis()->SetTitle("counts");
127  for (int i = 0; i < nPXDChips; i++) {
128  VxdID id = gTools->getChipIDFromPXDIndex(i);
129  int iLayer = id.getLayerNumber();
130  int iLadder = id.getLadderNumber();
131  int iSensor = id.getSensorNumber();
132  int iChip = gTools->getPXDChipNumber(id);
133  int IsU = gTools->isPXDSideU(id);
134  TString AxisTicks = Form("%i_%i_%i_u%iDCD", iLayer, iLadder, iSensor, iChip);
135  if (!IsU)
136  AxisTicks = Form("%i_%i_%i_v%iSWB", iLayer, iLadder, iSensor, iChip);
137  m_hitMapCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
138  m_hitMapClCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
139  }
140 
141  for (int i = 0; i < nPXDSensors; i++) {
142  VxdID id = gTools->getSensorIDFromPXDIndex(i);
143  int iLayer = id.getLayerNumber();
144  int iLadder = id.getLadderNumber();
145  int iSensor = id.getSensorNumber();
146  TString AxisTicks = Form("%i_%i_%i", iLayer, iLadder, iSensor);
147  m_hitMapCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
148  m_hitMapClCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
149  }
150 
151  m_fired.resize(nPXDSensors);
152  m_goodfired.resize(nPXDSensors);
153  m_clusters.resize(nPXDSensors);
154  m_goodclusters.resize(nPXDSensors);
155  m_startRow.resize(nPXDSensors);
156  m_chargStartRow.resize(nPXDSensors);
157  m_startRowCount.resize(nPXDSensors);
158  m_clusterCharge.resize(nPXDSensors);
159  m_clusterEnergy.resize(nPXDSensors);
160  m_pixelSignal.resize(nPXDSensors);
161  m_clusterSizeU.resize(nPXDSensors);
162  m_clusterSizeV.resize(nPXDSensors);
163  m_clusterSizeUV.resize(nPXDSensors);
164 
165  m_hitMapU.resize(nPXDSensors);
166  m_hitMapV.resize(nPXDSensors);
167  m_hitMap.resize(nPXDSensors);
168  m_hitMapUCl.resize(nPXDSensors);
169  m_hitMapVCl.resize(nPXDSensors);
170  m_hitMapCl.resize(nPXDSensors);
171  m_seed.resize(nPXDSensors);
172  for (int i = 0; i < nPXDSensors; i++) {
173  VxdID id = gTools->getSensorIDFromPXDIndex(i);
174  int iLayer = id.getLayerNumber();
175  int iLadder = id.getLadderNumber();
176  int iSensor = id.getSensorNumber();
177  VxdID sensorID(iLayer, iLadder, iSensor);
178  PXD::SensorInfo SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
179  string sensorDescr = str(format("%1%_%2%_%3%") % iLayer % iLadder % iSensor);
180  auto nUPixels = SensorInfo.getUCells();
181  auto nVPixels = SensorInfo.getVCells();
182  //----------------------------------------------------------------
183  // Number of fired pixels per frame
184  //----------------------------------------------------------------
185  string name = str(format("DQM_PXD_%1%_Fired") % sensorDescr);
186  string title = str(format("PXD Sensor %1% Fired pixels") % sensorDescr);
187  m_fired[i] = new TH1D(name.c_str(), title.c_str(), 50, 0, 50);
188  m_fired[i]->GetXaxis()->SetTitle("# of fired pixels");
189  m_fired[i]->GetYaxis()->SetTitle("counts");
190  //----------------------------------------------------------------
191  // Number of good fired pixels per frame
192  //----------------------------------------------------------------
193  name = str(format("DQM_PXD_%1%_GoodFired") % sensorDescr);
194  title = str(format("PXD Sensor %1% Good Fired pixels") % sensorDescr);
195  m_goodfired[i] = new TH1D(name.c_str(), title.c_str(), 50, 0, 50);
196  m_goodfired[i]->GetXaxis()->SetTitle("# of fired pixels");
197  m_goodfired[i]->GetYaxis()->SetTitle("counts");
198  //----------------------------------------------------------------
199  // Number of clusters per frame
200  //----------------------------------------------------------------
201  name = str(format("DQM_PXD_%1%_Clusters") % sensorDescr);
202  title = str(format("PXD Sensor %1% Number of clusters") % sensorDescr);
203  m_clusters[i] = new TH1D(name.c_str(), title.c_str(), 20, 0, 20);
204  m_clusters[i]->GetXaxis()->SetTitle("# of clusters");
205  m_clusters[i]->GetYaxis()->SetTitle("counts");
206  //----------------------------------------------------------------
207  // Number of good clusters per frame
208  //----------------------------------------------------------------
209  name = str(format("DQM_PXD_%1%_GoodClusters") % sensorDescr);
210  title = str(format("PXD Sensor %1% Number of good clusters") % sensorDescr);
211  m_goodclusters[i] = new TH1D(name.c_str(), title.c_str(), 20, 0, 20);
212  m_goodclusters[i]->GetXaxis()->SetTitle("# of clusters");
213  m_goodclusters[i]->GetYaxis()->SetTitle("counts");
214  //----------------------------------------------------------------
215  // Start row distribution
216  //----------------------------------------------------------------
217  name = str(format("DQM_PXD_%1%_StartRow") % sensorDescr);
218  title = str(format("PXD Sensor %1% Start row distribution") % sensorDescr);
219 
220  m_startRow[i] = new TH1D(name.c_str(), title.c_str(), nVPixels / 4, 0.0, nVPixels);
221  m_startRow[i]->GetXaxis()->SetTitle("start row [pitch units]");
222  m_startRow[i]->GetYaxis()->SetTitle("count");
223  //----------------------------------------------------------------
224  // Cluster seed charge by distance from the start row
225  //----------------------------------------------------------------
226  name = str(format("DQM_PXD_%1%_AverageSeedByStartRow") % sensorDescr);
227  title = str(format("PXD Sensor %1% Average seed charge by distance from the start row") % sensorDescr);
228  m_chargStartRow[i] = new TH1D(name.c_str(), title.c_str(), nVPixels / 4, 0.0, nVPixels);
229  m_chargStartRow[i]->GetXaxis()->SetTitle("distance from the start row [pitch units]");
230  m_chargStartRow[i]->GetYaxis()->SetTitle("average seed [ADU]");
231  name = str(format("DQM_PXD_%1%_SeedCountsByStartRow") % sensorDescr);
232  title = str(format("PXD Sensor %1% Seed charge count by distance from the start row") % sensorDescr);
233  m_startRowCount[i] = new TH1D(name.c_str(), title.c_str(), nVPixels / 4, 0.0, nVPixels);
234  m_startRowCount[i]->GetXaxis()->SetTitle("distance from the start row [pitch units]");
235  m_startRowCount[i]->GetYaxis()->SetTitle("count");
236  //----------------------------------------------------------------
237  // Cluster Charge
238  //----------------------------------------------------------------
239  name = str(format("DQM_PXD_%1%_ClusterCharge") % sensorDescr);
240  title = str(format("PXD Sensor %1% Cluster Charge") % sensorDescr);
241  m_clusterCharge[i] = new TH1D(name.c_str(), title.c_str(), 256, 0, 256);
242  m_clusterCharge[i]->GetXaxis()->SetTitle("charge of clusters [ADU]");
243  m_clusterCharge[i]->GetYaxis()->SetTitle("counts");
244  //----------------------------------------------------------------
245  // Cluster Energy
246  //----------------------------------------------------------------
247  name = str(format("DQM_PXD_%1%_ClusterEnergy") % sensorDescr);
248  title = str(format("PXD Sensor %1% Cluster Energy") % sensorDescr);
249  m_clusterEnergy[i] = new TH1D(name.c_str(), title.c_str(), 100, 0, 50);
250  m_clusterEnergy[i]->GetXaxis()->SetTitle("energy of clusters [keV]");
251  m_clusterEnergy[i]->GetYaxis()->SetTitle("counts");
252  //----------------------------------------------------------------
253  // Pixel Signal
254  //----------------------------------------------------------------
255  name = str(format("DQM_PXD_%1%_PixelSignal") % sensorDescr);
256  title = str(format("PXD Sensor %1% Pixel Signal") % sensorDescr);
257  m_pixelSignal[i] = new TH1D(name.c_str(), title.c_str(), 256, 0, 256);
258  m_pixelSignal[i]->GetXaxis()->SetTitle("signal of pixels [ADU]");
259  m_pixelSignal[i]->GetYaxis()->SetTitle("counts");
260  //----------------------------------------------------------------
261  // Cluster Size in U
262  //----------------------------------------------------------------
263  name = str(format("DQM_PXD_%1%_ClusterSizeU") % sensorDescr);
264  title = str(format("PXD Sensor %1% Cluster Size U") % sensorDescr);
265  m_clusterSizeU[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
266  m_clusterSizeU[i]->GetXaxis()->SetTitle("size of u clusters");
267  m_clusterSizeU[i]->GetYaxis()->SetTitle("counts");
268  //----------------------------------------------------------------
269  // Cluster Size in V
270  //----------------------------------------------------------------
271  name = str(format("DQM_PXD_%1%_ClusterSizeV") % sensorDescr);
272  title = str(format("PXD Sensor %1% Cluster Size V") % sensorDescr);
273  m_clusterSizeV[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
274  m_clusterSizeV[i]->GetXaxis()->SetTitle("size of v clusters");
275  m_clusterSizeV[i]->GetYaxis()->SetTitle("counts");
276  //----------------------------------------------------------------
277  // Cluster Size in U+V
278  //----------------------------------------------------------------
279  name = str(format("DQM_PXD_%1%_ClusterSizeUV") % sensorDescr);
280  title = str(format("PXD Sensor %1% Cluster Size U+V") % sensorDescr);
281  m_clusterSizeUV[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
282  m_clusterSizeUV[i]->GetXaxis()->SetTitle("size of u+v clusters");
283  m_clusterSizeUV[i]->GetYaxis()->SetTitle("counts");
284 
285  //----------------------------------------------------------------
286  // Hitmaps: Number of pixels by coordinate
287  //----------------------------------------------------------------
288  // Hitmaps in U
289  name = str(format("PXD_%1%_PixelHitmapU") % sensorDescr);
290  title = str(format("PXD Sensor %1% Pixel Hitmap in U") % sensorDescr);
291  m_hitMapU[i] = new TH1D(name.c_str(), title.c_str(), nUPixels, 0, nUPixels);
292  m_hitMapU[i]->GetXaxis()->SetTitle("u position [pitch units]");
293  m_hitMapU[i]->GetYaxis()->SetTitle("hits");
294  // Hitmaps in V
295  name = str(format("PXD_%1%_PixelHitmapV") % sensorDescr);
296  title = str(format("PXD Sensor %1% Pixel Hitmap in V") % sensorDescr);
297  m_hitMapV[i] = new TH1D(name.c_str(), title.c_str(), nVPixels, 0, nVPixels);
298  m_hitMapV[i]->GetXaxis()->SetTitle("v position [pitch units]");
299  m_hitMapV[i]->GetYaxis()->SetTitle("hits");
300  // Hitmaps in UV
301  name = str(format("PXD_%1%_PixelHitmap") % sensorDescr);
302  title = str(format("PXD Sensor %1% Pixel Hitmap") % sensorDescr);
303  m_hitMap[i] = new TH2D(name.c_str(), title.c_str(), nUPixels, 0, nUPixels, nVPixels, 0, nVPixels);
304  m_hitMap[i]->GetXaxis()->SetTitle("u position [pitch units]");
305  m_hitMap[i]->GetYaxis()->SetTitle("v position [pitch units]");
306  m_hitMap[i]->GetZaxis()->SetTitle("hits");
307 
308  //----------------------------------------------------------------
309  // Hitmaps: Number of clusters by coordinate
310  //----------------------------------------------------------------
311  // Hitmaps in U
312  name = str(format("PXD_%1%_HitmapClstU") % sensorDescr);
313  title = str(format("PXD Sensor %1% Hitmap Clusters in U") % sensorDescr);
314  m_hitMapUCl[i] = new TH1D(name.c_str(), title.c_str(), nUPixels, 0, nUPixels);
315  m_hitMapUCl[i]->GetXaxis()->SetTitle("u position [pitch units]");
316  m_hitMapUCl[i]->GetYaxis()->SetTitle("hits");
317  // Hitmaps in V
318  name = str(format("PXD_%1%_HitmapClstV") % sensorDescr);
319  title = str(format("PXD Sensor %1% Hitmap Clusters in V") % sensorDescr);
320  m_hitMapVCl[i] = new TH1D(name.c_str(), title.c_str(), nVPixels, 0, nVPixels);
321  m_hitMapVCl[i]->GetXaxis()->SetTitle("v position [pitch units]");
322  m_hitMapVCl[i]->GetYaxis()->SetTitle("hits");
323  // Hitmaps in UV
324  name = str(format("PXD_%1%_HitmapClst") % sensorDescr);
325  title = str(format("PXD Sensor %1% Hitmap Clusters") % sensorDescr);
326  m_hitMapCl[i] = new TH2D(name.c_str(), title.c_str(), nUPixels, 0, nUPixels, nVPixels, 0, nVPixels);
327  m_hitMapCl[i]->GetXaxis()->SetTitle("u position [pitch units]");
328  m_hitMapCl[i]->GetYaxis()->SetTitle("v position [pitch units]");
329  m_hitMapCl[i]->GetZaxis()->SetTitle("hits");
330 
331  //----------------------------------------------------------------
332  // Cluster seed charge distribution
333  //----------------------------------------------------------------
334  name = str(format("PXD_%1%_Seed") % sensorDescr);
335  title = str(format("PXD Sensor %1% Seed charge") % sensorDescr);
336  m_seed[i] = new TH1D(name.c_str(), title.c_str(), 256, 0, 256);
337  m_seed[i]->GetXaxis()->SetTitle("seed charge of clusters [ADU]");
338  m_seed[i]->GetYaxis()->SetTitle("count");
339 
340  }
341  oldDir->cd();
342 }
343 
344 
345 void PXDDQMClustersModule::initialize()
346 {
347  // Register histograms (calls back defineHisto)
348  REG_HISTOGRAM
349 
350  m_storeDAQEvtStats.isOptional();
351 
352  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
353  if (gTools->getNumberOfPXDLayers() != 0) {
354  //Register collections
355  StoreArray<PXDDigit> storePXDDigits(m_storePXDDigitsName);
356  StoreArray<PXDCluster> storePXDClusters(m_storePXDClustersName);
357  RelationArray relPXDClusterDigits(storePXDClusters, storePXDDigits);
358  m_storePXDClustersName = storePXDClusters.getName();
359  m_relPXDClusterDigitName = relPXDClusterDigits.getName();
360 
361  //Store names to speed up creation later
362  m_storePXDDigitsName = storePXDDigits.getName();
363  }
364 }
365 
366 void PXDDQMClustersModule::beginRun()
367 {
368  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
369  if (gTools->getNumberOfPXDLayers() == 0) return;
370 
371  if (m_hitMapCounts != nullptr) m_hitMapCounts->Reset();
372  if (m_hitMapFilterCounts != nullptr) m_hitMapFilterCounts->Reset();
373  if (m_hitMapClCounts != nullptr) m_hitMapClCounts->Reset();
374  if (m_hitMapClFilterCounts != nullptr) m_hitMapClFilterCounts->Reset();
375  if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Reset();
376  if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Reset();
377 
378  for (int i = 0; i < gTools->getNumberOfPXDSensors(); i++) {
379  if (m_fired[i] != nullptr) m_fired[i]->Reset();
380  if (m_goodfired[i] != nullptr) m_goodfired[i]->Reset();
381  if (m_clusters[i] != nullptr) m_clusters[i]->Reset();
382  if (m_goodclusters[i] != nullptr) m_goodclusters[i]->Reset();
383  if (m_startRow[i] != nullptr) m_startRow[i]->Reset();
384  if (m_chargStartRow[i] != nullptr) m_chargStartRow[i]->Reset();
385  if (m_startRowCount[i] != nullptr) m_startRowCount[i]->Reset();
386  if (m_clusterCharge[i] != nullptr) m_clusterCharge[i]->Reset();
387  if (m_clusterEnergy[i] != nullptr) m_clusterEnergy[i]->Reset();
388  if (m_pixelSignal[i] != nullptr) m_pixelSignal[i]->Reset();
389  if (m_clusterSizeU[i] != nullptr) m_clusterSizeU[i]->Reset();
390  if (m_clusterSizeV[i] != nullptr) m_clusterSizeV[i]->Reset();
391  if (m_clusterSizeUV[i] != nullptr) m_clusterSizeUV[i]->Reset();
392 
393  if (m_hitMapU[i] != nullptr) m_hitMapU[i]->Reset();
394  if (m_hitMapV[i] != nullptr) m_hitMapV[i]->Reset();
395  if (m_hitMap[i] != nullptr) m_hitMap[i]->Reset();
396  if (m_hitMapUCl[i] != nullptr) m_hitMapUCl[i]->Reset();
397  if (m_hitMapVCl[i] != nullptr) m_hitMapVCl[i]->Reset();
398  if (m_hitMapCl[i] != nullptr) m_hitMapCl[i]->Reset();
399  if (m_seed[i] != nullptr) m_seed[i]->Reset();
400 
401  }
402 }
403 
404 
405 void PXDDQMClustersModule::event()
406 {
407  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
408  if (gTools->getNumberOfPXDLayers() == 0) return;
409 
410  const StoreArray<PXDDigit> storePXDDigits(m_storePXDDigitsName);
411  const StoreArray<PXDCluster> storePXDClusters(m_storePXDClustersName);
412  const RelationArray relPXDClusterDigits(storePXDClusters, storePXDDigits, m_relPXDClusterDigitName);
413 
414  // If there are no digits, leave
415  if (!storePXDDigits || !storePXDDigits.getEntries()) return;
416 
417  int firstPXDLayer = gTools->getFirstPXDLayer();
418  int lastPXDLayer = gTools->getLastPXDLayer();
419  int nPXDSensors = gTools->getNumberOfPXDSensors();
420 
421  // PXD basic histograms:
422  // Fired pixels
423  vector< int > Pixels(nPXDSensors);
424  vector< int > GoodPixels(nPXDSensors);
425  for (const PXDDigit& digit : storePXDDigits) {
426  int iLayer = digit.getSensorID().getLayerNumber();
427  if ((iLayer < firstPXDLayer) || (iLayer > lastPXDLayer)) continue;
428  int iLadder = digit.getSensorID().getLadderNumber();
429  int iSensor = digit.getSensorID().getSensorNumber();
430  VxdID sensorID(iLayer, iLadder, iSensor);
431  int index = gTools->getPXDSensorIndex(sensorID);
432  PXD::SensorInfo SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
433  Pixels[index]++;
434  int iChip = PXDMappingLookup::getDCDID(digit.getUCellID(), digit.getVCellID(), sensorID);
435  int indexChip = gTools->getPXDChipIndex(sensorID, kTRUE, iChip);
436  if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Fill(indexChip);
437  iChip = PXDMappingLookup::getSWBID(digit.getVCellID());
438  indexChip = gTools->getPXDChipIndex(sensorID, kFALSE, iChip);
439  if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Fill(indexChip);
440  if (m_pixelSignal[index] != nullptr) m_pixelSignal[index]->Fill(digit.getCharge());
441  if (m_hitMapCounts != nullptr) m_hitMapCounts->Fill(index);
442  // filter for pixel charge
443  if (digit.getCharge() >= m_CutMinCharge && digit.getCharge() < 255) {
444  GoodPixels[index]++;
445  if (m_hitMapFilterCounts != nullptr) m_hitMapFilterCounts->Fill(index);
446  if (m_hitMapU[index] != nullptr) m_hitMapU[index]->Fill(digit.getUCellID());
447  if (m_hitMapV[index] != nullptr) m_hitMapV[index]->Fill(digit.getVCellID());
448  if (m_hitMap[index] != nullptr) m_hitMap[index]->Fill(digit.getUCellID(), digit.getVCellID());
449  }
450  }
451  for (int i = 0; i < nPXDSensors; i++) {
452  if (m_fired[i] != nullptr && Pixels[i] > 0) m_fired[i]->Fill(Pixels[i]);
453  if (m_goodfired[i] != nullptr && GoodPixels[i] > 0) m_goodfired[i]->Fill(GoodPixels[i]);
454  }
455 
456  // Hitmaps, Charge, Seed, Size, ...
457  vector< int > Clusters(nPXDSensors);
458  vector< int > GoodClusters(nPXDSensors);
459  for (const PXDCluster& cluster : storePXDClusters) {
460  int iLayer = cluster.getSensorID().getLayerNumber();
461  if ((iLayer < firstPXDLayer) || (iLayer > lastPXDLayer)) continue;
462  int iLadder = cluster.getSensorID().getLadderNumber();
463  int iSensor = cluster.getSensorID().getSensorNumber();
464  VxdID sensorID(iLayer, iLadder, iSensor);
465  int index = gTools->getPXDSensorIndex(sensorID);
466  PXD::SensorInfo SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
467  Clusters[index]++;
468  int iChip = PXDMappingLookup::getDCDID(SensorInfo.getUCellID(cluster.getU()), SensorInfo.getVCellID(cluster.getV()), sensorID);
469  int indexChip = gTools->getPXDChipIndex(sensorID, kTRUE, iChip);
470  if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Fill(indexChip);
471  iChip = PXDMappingLookup::getSWBID(SensorInfo.getVCellID(cluster.getV()));
472  indexChip = gTools->getPXDChipIndex(sensorID, kFALSE, iChip);
473  if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Fill(indexChip);
474  if (m_hitMapClCounts != nullptr) m_hitMapClCounts->Fill(index);
475  if (m_clusterCharge[index] != nullptr) m_clusterCharge[index]->Fill(cluster.getCharge());
476  // FIXME: The cluster charge is stored in ADU. We have to extract the
477  // area dependent conversion factor ADU->eV. Assume this is the same for all pixels of the
478  // cluster.
479  auto cluster_uID = SensorInfo.getUCellID(cluster.getU());
480  auto cluster_vID = SensorInfo.getVCellID(cluster.getV());
481  auto ADUToEnergy = PXDGainCalibrator::getInstance().getADUToEnergy(sensorID, cluster_uID, cluster_vID);
482  if (m_clusterEnergy[index] != nullptr) m_clusterEnergy[index]->Fill(cluster.getCharge()* ADUToEnergy / Unit::keV);
483  if (m_clusterSizeU[index] != nullptr) m_clusterSizeU[index]->Fill(cluster.getUSize());
484  if (m_clusterSizeV[index] != nullptr) m_clusterSizeV[index]->Fill(cluster.getVSize());
485  if (m_clusterSizeUV[index] != nullptr) m_clusterSizeUV[index]->Fill(cluster.getSize());
486  if (m_seed[index] != nullptr) m_seed[index]->Fill(cluster.getSeedCharge());
487 
488  // filter for cluster size and seed pixel
489  if (cluster.getSize() <= m_CutMaxClusterSize && cluster.getCharge() >= m_CutMinClusterCharge
490  && cluster.getSeedCharge() >= m_CutMinSeedCharge && cluster.getSeedCharge() < 255) {
491  GoodClusters[index]++;
492  if (m_hitMapClFilterCounts != nullptr) m_hitMapClFilterCounts->Fill(index);
493 
494  if (m_hitMapUCl[index] != nullptr) m_hitMapUCl[index]->Fill(
495  SensorInfo.getUCellID(cluster.getU()));
496  if (m_hitMapVCl[index] != nullptr) m_hitMapVCl[index]->Fill(
497  SensorInfo.getVCellID(cluster.getV()));
498  if (m_hitMapCl[index] != nullptr) m_hitMapCl[index]->Fill(
499  SensorInfo.getUCellID(cluster.getU()),
500  SensorInfo.getVCellID(cluster.getV()));
501  }
502  }
503  for (int i = 0; i < nPXDSensors; i++) {
504  if (m_clusters[i] != nullptr && Clusters[i] > 0) m_clusters[i]->Fill(Clusters[i]);
505  if (m_goodclusters[i] != nullptr && GoodClusters[i] > 0) m_goodclusters[i]->Fill(GoodClusters[i]);
506  }
507 
508  // Only fill start row (triggergate) histos when data is available
509  if (m_storeDAQEvtStats.isValid()) {
510  std::map<VxdID, unsigned short> startRows;
511  for (int index = 0; index < nPXDSensors; index++) {
512  VxdID id = gTools->getSensorIDFromPXDIndex(index);
513 
514  const PXDDAQDHEStatus* dhe = (*m_storeDAQEvtStats).findDHE(id);
515  if (dhe != nullptr) {
516  auto startRow = dhe->getStartRow();
517  if (m_startRow[index] != nullptr) m_startRow[index]->Fill(startRow);
518  startRows.insert(std::make_pair(id, startRow));
519  } else {
520  B2WARNING("No PXDDAQDHEStatus for VXD Sensor " << id << " found.");
521  }
522  }
523 
524  // Cluster seed charge by start row
525  for (auto& cluster : storePXDClusters) {
526  VxdID sensorID = cluster.getSensorID();
527  int index = gTools->getPXDSensorIndex(sensorID);
528  PXD::SensorInfo SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
529 
530  float fDistance = SensorInfo.getVCellID(cluster.getV()) - startRows[cluster.getSensorID()];
531  if (fDistance < 0) fDistance += SensorInfo.getVCells();
532  if (m_chargStartRow[index] != nullptr) m_chargStartRow[index]->Fill(fDistance, cluster.getSeedCharge());
533  if (m_startRowCount[index] != nullptr) m_startRowCount[index]->Fill(fDistance);
534  }
535  }
536 
537 }
Belle2::VXD::SensorInfoBase::getUCells
int getUCells() const
Return number of pixel/strips in u direction.
Definition: SensorInfoBase.h:223
Belle2::RelationArray
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:72
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::StoreAccessorBase::getName
const std::string & getName() const
Return name under which the object is saved in the DataStore.
Definition: StoreAccessorBase.h:130
Belle2::PXDDQMClustersModule
PXD DQM Module.
Definition: PXDDQMClustersModule.h:41
Belle2::PXDDigit
The PXD digit class.
Definition: PXDDigit.h:38
Belle2::PXD::SensorInfo
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:34
Belle2::PXD
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Definition: PXDCalibrationUtilities.h:28
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::PXDDAQDHEStatus::getStartRow
unsigned short getStartRow(void) const
get Trigger Start Row
Definition: PXDDAQDHEStatus.h:122
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::PXDCluster
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:41
Belle2::PXDDAQDHEStatus
The PXD DAQ DHE Status class.
Definition: PXDDAQDHEStatus.h:46
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