Belle II Software development
PXDDQMExpressRecoModule.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/PXDDQMExpressRecoModule.h"
10
11#include <pxd/geometry/SensorInfo.h>
12#include <vxd/geometry/GeoCache.h>
13#include <vxd/geometry/SensorInfoBase.h>
14#include <vxd/geometry/GeoTools.h>
15#include <pxd/unpacking/PXDMappingLookup.h>
16
17#include <boost/format.hpp>
18
19#include "TDirectory.h"
20
21using namespace std;
22using boost::format;
23using namespace Belle2;
24using namespace Belle2::PXD;
25
26//-----------------------------------------------------------------
27// Register the Module
28//-----------------------------------------------------------------
29REG_MODULE(PXDDQMExpressReco);
30
31
32//-----------------------------------------------------------------
33// Implementation
34//-----------------------------------------------------------------
35
37{
38 //Set module properties
39 setDescription("PXD DQM module for Express Reco "
40 "Recommended Number of events for monitor is 40 kEvents or more to fill all histograms "
41 );
42
43 setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
44 addParam("CutMinCharge", m_CutMinCharge,
45 "cut on pixel charge for accepting good pixel, default >= 12", 12);
46 addParam("CutMinSeedCharge", m_CutMinSeedCharge,
47 "cut on cluster seed for accepting good cluster, default >= 12", 12);
48 addParam("CutMaxClusterSize", m_CutMaxClusterSize,
49 "cut on cluster size accepting good cluster, default <= 4", 4);
50 addParam("CutMinClusterCharge", m_CutMinClusterCharge,
51 "cut on cluster charge for accepting good cluster, default >= 12", 12);
52 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
53 std::string("PXDER"));
54}
55
56
57//------------------------------------------------------------------
58// Function to define histograms
59//-----------------------------------------------------------------
60
62{
64 if (gTools->getNumberOfLayers() == 0) {
65 B2FATAL("Missing geometry for VXD, check steering file.");
66 }
67 if (gTools->getNumberOfPXDLayers() == 0) {
68 B2WARNING("Missing geometry for PXD, PXD-DQM is skipped.");
69 return;
70 }
71
72 // Create a separate histogram directories and cd into it.
73 TDirectory* oldDir = gDirectory;
74 if (m_histogramDirectoryName != "") {
75 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
76 oldDir->cd(m_histogramDirectoryName.c_str());
77 }
78
79 // basic constants presets:
80 int nPXDSensors = gTools->getNumberOfPXDSensors();
81 int nPXDChips = gTools->getTotalPXDChips();
82
83 // Create basic histograms:
84 m_hitMapCounts = new TH1D("DQMER_PXD_PixelHitmapCounts", "PXD Integrated number of fired pixels per sensor",
85 nPXDSensors, 0, nPXDSensors);
86 m_hitMapCounts->GetXaxis()->SetTitle("Sensor ID");
87 m_hitMapCounts->GetYaxis()->SetTitle("counts");
88
89 m_hitMapFilterCounts = new TH1D("DQMER_PXD_PixelHitmapFilterCounts", "PXD Integrated number of filtered pixels per sensor",
90 nPXDSensors, 0, nPXDSensors);
91 m_hitMapFilterCounts->GetXaxis()->SetTitle("Sensor ID");
92 m_hitMapFilterCounts->GetYaxis()->SetTitle("counts");
93
94 m_hitMapClCounts = new TH1D("DQMER_PXD_ClusterHitmapCounts", "PXD Integrated number of clusters per sensor",
95 nPXDSensors, 0, nPXDSensors);
96 m_hitMapClCounts->GetXaxis()->SetTitle("Sensor ID");
97 m_hitMapClCounts->GetYaxis()->SetTitle("counts");
98
99 m_hitMapClFilterCounts = new TH1D("DQMER_PXD_ClusterHitmapFilterCounts", "PXD Integrated number of filtered clusters per sensor",
100 nPXDSensors, 0, nPXDSensors);
101 m_hitMapClFilterCounts->GetXaxis()->SetTitle("Sensor ID");
102 m_hitMapClFilterCounts->GetYaxis()->SetTitle("counts");
103
104 // basic counters per chip:
105 m_hitMapCountsChip = new TH1D("DQMER_PXD_PixelHitmapCountsChip", "PXD Integrated number of fired pixels per chip",
106 nPXDChips, 0, nPXDChips);
107 m_hitMapCountsChip->GetXaxis()->SetTitle("Chip ID");
108 m_hitMapCountsChip->GetYaxis()->SetTitle("counts");
109 m_hitMapClCountsChip = new TH1D("DQMER_PXD_ClusterHitmapCountsChip", "PXD Integrated number of clusters per chip",
110 nPXDChips, 0, nPXDChips);
111 m_hitMapClCountsChip->GetXaxis()->SetTitle("Chip ID");
112 m_hitMapClCountsChip->GetYaxis()->SetTitle("counts");
113 for (int i = 0; i < nPXDChips; i++) {
114 VxdID id = gTools->getChipIDFromPXDIndex(i);
115 int iLayer = id.getLayerNumber();
116 int iLadder = id.getLadderNumber();
117 int iSensor = id.getSensorNumber();
118 int iChip = gTools->getPXDChipNumber(id);
119 int IsU = gTools->isPXDSideU(id);
120 TString AxisTicks = Form("%i_%i_%i_u%iDCD", iLayer, iLadder, iSensor, iChip);
121 if (!IsU)
122 AxisTicks = Form("%i_%i_%i_v%iSWB", iLayer, iLadder, iSensor, iChip);
123 m_hitMapCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
124 m_hitMapClCountsChip->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
125 }
126
127 for (int i = 0; i < nPXDSensors; i++) {
128 VxdID id = gTools->getSensorIDFromPXDIndex(i);
129 int iLayer = id.getLayerNumber();
130 int iLadder = id.getLadderNumber();
131 int iSensor = id.getSensorNumber();
132 TString AxisTicks = Form("%i_%i_%i", iLayer, iLadder, iSensor);
133 m_hitMapCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
134 m_hitMapClCounts->GetXaxis()->SetBinLabel(i + 1, AxisTicks.Data());
135 }
136
137 m_fired.resize(nPXDSensors);
138 m_goodfired.resize(nPXDSensors);
139 m_clusters.resize(nPXDSensors);
140 m_goodclusters.resize(nPXDSensors);
141 // FIXME: startrow histos are for experts not shifters
142 //m_startRow.resize(nPXDSensors);
143 //m_chargStartRow.resize(nPXDSensors);
144 //m_startRowCount.resize(nPXDSensors);
145 m_clusterCharge.resize(nPXDSensors);
146 m_pixelSignal.resize(nPXDSensors);
147 m_clusterSizeU.resize(nPXDSensors);
148 m_clusterSizeV.resize(nPXDSensors);
149 m_clusterSizeUV.resize(nPXDSensors);
150 for (int i = 0; i < nPXDSensors; i++) {
151 VxdID id = gTools->getSensorIDFromPXDIndex(i);
152 int iLayer = id.getLayerNumber();
153 int iLadder = id.getLadderNumber();
154 int iSensor = id.getSensorNumber();
155 VxdID sensorID(iLayer, iLadder, iSensor);
157 string sensorDescr = str(format("%1%_%2%_%3%") % iLayer % iLadder % iSensor);
158 //----------------------------------------------------------------
159 // Number of fired pixels per frame
160 //----------------------------------------------------------------
161 string name = str(format("DQMER_PXD_%1%_Fired") % sensorDescr);
162 string title = str(format("PXD Sensor %1% Fired pixels") % sensorDescr);
163 m_fired[i] = new TH1D(name.c_str(), title.c_str(), 200, 0, 200);
164 m_fired[i]->SetCanExtend(TH1::kAllAxes);
165 m_fired[i]->GetXaxis()->SetTitle("# of fired pixels");
166 m_fired[i]->GetYaxis()->SetTitle("counts");
167 //----------------------------------------------------------------
168 // Number of good fired pixels per frame
169 //----------------------------------------------------------------
170 name = str(format("DQMER_PXD_%1%_GoodFired") % sensorDescr);
171 title = str(format("PXD Sensor %1% Good pixels") % sensorDescr);
172 m_goodfired[i] = new TH1D(name.c_str(), title.c_str(), 200, 0, 200);
173 m_goodfired[i]->SetCanExtend(TH1::kAllAxes);
174 m_goodfired[i]->GetXaxis()->SetTitle("# of fired pixels");
175 m_goodfired[i]->GetYaxis()->SetTitle("counts");
176 //----------------------------------------------------------------
177 // Number of clusters per frame
178 //----------------------------------------------------------------
179 name = str(format("DQMER_PXD_%1%_Clusters") % sensorDescr);
180 title = str(format("PXD Sensor %1% Clusters") % sensorDescr);
181 m_clusters[i] = new TH1D(name.c_str(), title.c_str(), 200, 0, 200);
182 m_clusters[i]->SetCanExtend(TH1::kAllAxes);
183 m_clusters[i]->GetXaxis()->SetTitle("# of clusters");
184 m_clusters[i]->GetYaxis()->SetTitle("counts");
185 //----------------------------------------------------------------
186 // Number of good clusters per frame
187 //----------------------------------------------------------------
188 name = str(format("DQMER_PXD_%1%_GoodClusters") % sensorDescr);
189 title = str(format("PXD Sensor %1% Good clusters") % sensorDescr);
190 m_goodclusters[i] = new TH1D(name.c_str(), title.c_str(), 200, 0, 200);
191 m_goodclusters[i]->SetCanExtend(TH1::kAllAxes);
192 m_goodclusters[i]->GetXaxis()->SetTitle("# of clusters");
193 m_goodclusters[i]->GetYaxis()->SetTitle("counts");
194 //----------------------------------------------------------------
195 // Start row distribution
196 // FIXME: expert level, remove here at some point
197 //----------------------------------------------------------------
198 //name = str(format("DQMER_PXD_%1%_StartRow") % sensorDescr);
199 //title = str(format("PXD Sensor %1% Start row distribution") % sensorDescr);
200
201 //int nPixels;/** Number of pixels on PXD v direction */
202 //nPixels = SensorInfo.getVCells();
203 //m_startRow[i] = new TH1D(name.c_str(), title.c_str(), nPixels / 4, 0.0, nPixels);
204 //m_startRow[i]->GetXaxis()->SetTitle("start row [pitch units]");
205 //m_startRow[i]->GetYaxis()->SetTitle("count");
206 //----------------------------------------------------------------
207 // Cluster seed charge by distance from the start row
208 //----------------------------------------------------------------
209 //name = str(format("DQMER_PXD_%1%_AverageSeedByStartRow") % sensorDescr);
210 //title = str(format("PXD Sensor %1% Average seed charge by distance from the start row") % sensorDescr);
211 //m_chargStartRow[i] = new TH1D(name.c_str(), title.c_str(), nPixels / 4, 0.0, nPixels);
212 //m_chargStartRow[i]->GetXaxis()->SetTitle("distance from the start row [pitch units]");
213 //m_chargStartRow[i]->GetYaxis()->SetTitle("average seed [ADU]");
214 //name = str(format("DQMER_PXD_%1%_SeedCountsByStartRow") % sensorDescr);
215 //title = str(format("PXD Sensor %1% Seed charge count by distance from the start row") % sensorDescr);
216 //m_startRowCount[i] = new TH1D(name.c_str(), title.c_str(), nPixels / 4, 0.0, nPixels);
217 //m_startRowCount[i]->GetXaxis()->SetTitle("distance from the start row [pitch units]");
218 //m_startRowCount[i]->GetYaxis()->SetTitle("count");
219 //----------------------------------------------------------------
220 // Cluster Charge
221 //----------------------------------------------------------------
222 name = str(format("DQMER_PXD_%1%_ClusterCharge") % sensorDescr);
223 title = str(format("PXD Sensor %1% Cluster Charge") % sensorDescr);
224 m_clusterCharge[i] = new TH1D(name.c_str(), title.c_str(), 256, 0, 256);
225 m_clusterCharge[i]->GetXaxis()->SetTitle("charge of clusters [ADU]");
226 m_clusterCharge[i]->GetYaxis()->SetTitle("counts");
227 //----------------------------------------------------------------
228 // Pixel Signal
229 //----------------------------------------------------------------
230 name = str(format("DQMER_PXD_%1%_PixelSignal") % sensorDescr);
231 title = str(format("PXD Sensor %1% Pixel Signal") % sensorDescr);
232 m_pixelSignal[i] = new TH1D(name.c_str(), title.c_str(), 256, 0, 256);
233 m_pixelSignal[i]->GetXaxis()->SetTitle("signal of pixels [ADU]");
234 m_pixelSignal[i]->GetYaxis()->SetTitle("counts");
235 //----------------------------------------------------------------
236 // Cluster Size in U
237 //----------------------------------------------------------------
238 name = str(format("DQMER_PXD_%1%_ClusterSizeU") % sensorDescr);
239 title = str(format("PXD Sensor %1% Cluster Size U") % sensorDescr);
240 m_clusterSizeU[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
241 m_clusterSizeU[i]->GetXaxis()->SetTitle("size of u clusters");
242 m_clusterSizeU[i]->GetYaxis()->SetTitle("counts");
243 //----------------------------------------------------------------
244 // Cluster Size in V
245 //----------------------------------------------------------------
246 name = str(format("DQMER_PXD_%1%_ClusterSizeV") % sensorDescr);
247 title = str(format("PXD Sensor %1% Cluster Size V") % sensorDescr);
248 m_clusterSizeV[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
249 m_clusterSizeV[i]->GetXaxis()->SetTitle("size of v clusters");
250 m_clusterSizeV[i]->GetYaxis()->SetTitle("counts");
251 //----------------------------------------------------------------
252 // Cluster Size in U+V
253 //----------------------------------------------------------------
254 name = str(format("DQMER_PXD_%1%_ClusterSizeUV") % sensorDescr);
255 title = str(format("PXD Sensor %1% Cluster Size U+V") % sensorDescr);
256 m_clusterSizeUV[i] = new TH1D(name.c_str(), title.c_str(), 10, 0, 10);
257 m_clusterSizeUV[i]->GetXaxis()->SetTitle("size of u+v clusters");
258 m_clusterSizeUV[i]->GetYaxis()->SetTitle("counts");
259 }
260
261 oldDir->cd();
262}
263
264
266{
267 // Register histograms (calls back defineHisto)
268 REG_HISTOGRAM
269
270 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
271 if (gTools->getNumberOfPXDLayers() != 0) {
272 //Register collections
275 }
276}
277
279{
280 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
281 if (gTools->getNumberOfPXDLayers() == 0) return;
282
283 if (m_hitMapCounts != nullptr) m_hitMapCounts->Reset();
284 if (m_hitMapFilterCounts != nullptr) m_hitMapFilterCounts->Reset();
285 if (m_hitMapClCounts != nullptr) m_hitMapClCounts->Reset();
286 if (m_hitMapClFilterCounts != nullptr) m_hitMapClFilterCounts->Reset();
287 if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Reset();
288 if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Reset();
289
290 for (int i = 0; i < gTools->getNumberOfPXDSensors(); i++) {
291 if (m_fired[i] != nullptr) m_fired[i]->Reset();
292 if (m_goodfired[i] != nullptr) m_goodfired[i]->Reset();
293 if (m_clusters[i] != nullptr) m_clusters[i]->Reset();
294 if (m_goodclusters[i] != nullptr) m_goodclusters[i]->Reset();
295 // FIXME: startrow is for expert only, not shifter
296 //if (m_startRow[i] != nullptr) m_startRow[i]->Reset();
297 //if (m_chargStartRow[i] != nullptr) m_chargStartRow[i]->Reset();
298 //if (m_startRowCount[i] != nullptr) m_startRowCount[i]->Reset();
299 if (m_clusterCharge[i] != nullptr) m_clusterCharge[i]->Reset();
300 if (m_pixelSignal[i] != nullptr) m_pixelSignal[i]->Reset();
301 if (m_clusterSizeU[i] != nullptr) m_clusterSizeU[i]->Reset();
302 if (m_clusterSizeV[i] != nullptr) m_clusterSizeV[i]->Reset();
303 if (m_clusterSizeUV[i] != nullptr) m_clusterSizeUV[i]->Reset();
304 }
305
306}
307
308
310{
311 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
312 if (gTools->getNumberOfPXDLayers() == 0) return;
313
314 // If there are no digits, leave
315 if (!m_storePXDDigits || !m_storePXDDigits.getEntries()) return;
316
317 int nPXDSensors = gTools->getNumberOfPXDSensors();
318
319 // PXD basic histograms:
320 // Fired pixels
321 vector< int > Pixels(nPXDSensors);
322 vector< int > GoodPixels(nPXDSensors);
323 for (const PXDDigit& digit : m_storePXDDigits) {
324 int iLayer = digit.getSensorID().getLayerNumber();
325 int iLadder = digit.getSensorID().getLadderNumber();
326 int iSensor = digit.getSensorID().getSensorNumber();
327 VxdID sensorID(iLayer, iLadder, iSensor);
328 int index = gTools->getPXDSensorIndex(sensorID);
330 Pixels[index]++;
331 int iChip = PXDMappingLookup::getDCDID(digit.getUCellID(), digit.getVCellID(), sensorID);
332 int indexChip = gTools->getPXDChipIndex(sensorID, kTRUE, iChip);
333 if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Fill(indexChip);
334 iChip = PXDMappingLookup::getSWBID(digit.getVCellID());
335 indexChip = gTools->getPXDChipIndex(sensorID, kFALSE, iChip);
336 if (m_hitMapCountsChip != nullptr) m_hitMapCountsChip->Fill(indexChip);
337 if (m_pixelSignal[index] != nullptr) m_pixelSignal[index]->Fill(digit.getCharge());
338 if (m_hitMapCounts != nullptr) m_hitMapCounts->Fill(index);
339 // filter for pixel charge
340 if (digit.getCharge() >= m_CutMinCharge && digit.getCharge() < 255) {
341 GoodPixels[index]++;
342 if (m_hitMapFilterCounts != nullptr) m_hitMapFilterCounts->Fill(index);
343 }
344 }
345 for (int i = 0; i < nPXDSensors; i++) {
346 if (m_fired[i] != nullptr && Pixels[i] > 0) m_fired[i]->Fill(Pixels[i]);
347 if (m_goodfired[i] != nullptr && GoodPixels[i] > 0) m_goodfired[i]->Fill(GoodPixels[i]);
348 }
349
350 // Hitmaps, Charge, Seed, Size, ...
351 vector< int > Clusters(nPXDSensors);
352 vector< int > GoodClusters(nPXDSensors);
353 for (const PXDCluster& cluster : m_storePXDClusters) {
354 int iLayer = cluster.getSensorID().getLayerNumber();
355 int iLadder = cluster.getSensorID().getLadderNumber();
356 int iSensor = cluster.getSensorID().getSensorNumber();
357 VxdID sensorID(iLayer, iLadder, iSensor);
358 int index = gTools->getPXDSensorIndex(sensorID);
360 Clusters[index]++;
361 int iChip = PXDMappingLookup::getDCDID(SensorInfo.getUCellID(cluster.getU()), SensorInfo.getVCellID(cluster.getV()), sensorID);
362 int indexChip = gTools->getPXDChipIndex(sensorID, kTRUE, iChip);
363 if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Fill(indexChip);
364 iChip = PXDMappingLookup::getSWBID(SensorInfo.getVCellID(cluster.getV()));
365 indexChip = gTools->getPXDChipIndex(sensorID, kFALSE, iChip);
366 if (m_hitMapClCountsChip != nullptr) m_hitMapClCountsChip->Fill(indexChip);
367 if (m_hitMapClCounts != nullptr) m_hitMapClCounts->Fill(index);
368 if (m_clusterCharge[index] != nullptr) m_clusterCharge[index]->Fill(cluster.getCharge());
369 if (m_clusterSizeU[index] != nullptr) m_clusterSizeU[index]->Fill(cluster.getUSize());
370 if (m_clusterSizeV[index] != nullptr) m_clusterSizeV[index]->Fill(cluster.getVSize());
371 if (m_clusterSizeUV[index] != nullptr) m_clusterSizeUV[index]->Fill(cluster.getSize());
372 // filter for cluster size and seed pixel
373 if (cluster.getSize() <= m_CutMaxClusterSize && cluster.getCharge() >= m_CutMinClusterCharge
374 && cluster.getSeedCharge() >= m_CutMinSeedCharge && cluster.getSeedCharge() < 255) {
375 GoodClusters[index]++;
376 if (m_hitMapClFilterCounts != nullptr) m_hitMapClFilterCounts->Fill(index);
377 }
378 }
379 for (int i = 0; i < nPXDSensors; i++) {
380 if (m_clusters[i] != nullptr && Clusters[i] > 0) m_clusters[i]->Fill(Clusters[i]);
381 if (m_goodclusters[i] != nullptr && GoodClusters[i] > 0) m_goodclusters[i]->Fill(GoodClusters[i]);
382 }
383}
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
std::vector< TH1D * > m_goodfired
Filtered fired pixels per event.
TH1D * m_hitMapClFilterCounts
Hitmaps of filtered Clusters.
void initialize() override final
Initialize.
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.
std::string m_storePXDClustersName
PXDClusters StoreArray name.
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_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_goodclusters
filtered Clusters per event
StoreArray< PXDDigit > m_storePXDDigits
Storearray for Digits
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
StoreArray< PXDCluster > m_storePXDClusters
Storearray for Cluster
TH1D * m_hitMapClCountsChip
Hitmaps of clusters on chips.
std::vector< TH1D * > m_clusterCharge
Start row distribution.
std::string m_storePXDDigitsName
PXDDigits StoreArray name.
The PXD digit class.
Definition: PXDDigit.h:27
static int getDCDID(const int u, const int v, const VxdID sensorID)
get ID of DCD for giving pixel, range: 1..4.
static int getSWBID(const int v)
get ID of SWB for giving pixel, range: 1..6.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:142
unsigned short getNumberOfPXDSensors() const
Get number of PXD sensors.
Definition: GeoTools.h:133
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.
STL namespace.