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