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