Belle II Software  release-05-02-19
PXDHotPixelMaskCalibrationAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Benjamin Schwenker *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <pxd/calibration/PXDHotPixelMaskCalibrationAlgorithm.h>
12 #include <pxd/dbobjects/PXDMaskedPixelPar.h>
13 #include <pxd/dbobjects/PXDDeadPixelPar.h>
14 #include <pxd/dbobjects/PXDOccupancyInfoPar.h>
15 
16 #include <boost/format.hpp>
17 #include <string>
18 #include <vector>
19 #include "TH1I.h"
20 #include "TMath.h"
21 
22 using namespace std;
23 using boost::format;
24 using namespace Belle2;
25 
26 
27 PXDHotPixelMaskCalibrationAlgorithm::PXDHotPixelMaskCalibrationAlgorithm(): CalibrationAlgorithm("PXDHotPixelMaskCollector"),
28  forceContinueMasking(false), minEvents(10000), minHits(20), pixelMultiplier(7), maskDrains(true),
29  drainMultiplier(7), maskRows(true), rowMultiplier(7)
30 {
32  " -------------------------- PXDHotPixelMak Calibration Algorithm ------------------------\n"
33  " \n"
34  " Algorithm which masks hot pixels with too large occupancy and dead pixels w/o no hits. \n"
35  " ----------------------------------------------------------------------------------------\n"
36  );
37 }
38 
40 {
41  // Run list of current IoV and all runs to be calibrated
42  auto expRuns = getRunList();
43  auto expRunsAll = getRunListFromAllData();
44 
45  // Save the current directory to change back later
46  TDirectory* currentDir = gDirectory;
47 
48  // Create TFile if not exist
49  if (!m_file) {
50  std::string fileName = (this->getPrefix()) + "debug.root";
51  B2INFO("Creating file " << fileName);
52  m_file = std::make_shared<TFile>(fileName.c_str(), "RECREATE");
53  }
54 
55  //m_file->cd();
56  string iov_str = str(format("E%1%_R%2%_E%3%_R%4%") % expRuns.front().first
57  % expRuns.front().second % expRuns.back().first % expRuns.back().second);
58  m_file->mkdir(iov_str.c_str());
59  m_file->cd(iov_str.c_str());
60 
61  // Collector info
62  auto collector_pxdhits = getObjectPtr<TH1I>("PXDHits");
63  auto collector_pxdhitcounts = getObjectPtr<TH1I>("PXDHitCounts");
64 
65  auto nevents = collector_pxdhits->GetEntries();
66 
67  for (auto sensBin = 1; sensBin <= collector_pxdhitcounts->GetXaxis()->GetNbins(); sensBin++) {
68  string sensorDescr = collector_pxdhitcounts->GetXaxis()->GetBinLabel(sensBin);
69  VxdID id(sensorDescr);
70 
71  // Get hitmap from collector
72  string name = str(format("PXD_%1%_PixelHitmap") % id.getID());
73  auto collector_pxdhitmap = getObjectPtr<TH1I>(name.c_str());
74  if (collector_pxdhitmap == nullptr) {
75  continue;
76  }
77  TH1I* debugHisto = new TH1I(*collector_pxdhitmap);
78  // Set overflow to total number of events
79  debugHisto->SetBinContent(debugHisto->GetNbinsX() + 1, nevents);
80  debugHisto->Write();
81  delete debugHisto;
82  }
83 
84  // Close TFile
85  if (expRuns.back() == expRunsAll.back()) {
86  B2INFO("Reached Final ExpRun: (" << expRuns.back().first << ", " << expRuns.back().second << ")");
87  m_file->Close();
88 
89  }
90  currentDir->cd();
91 }
92 
94 {
95 
96  auto collector_pxdhits = getObjectPtr<TH1I>("PXDHits");
97  auto collector_pxdhitcounts = getObjectPtr<TH1I>("PXDHitCounts");
98 
99  // We should have some minimum number of events
100  auto nevents = collector_pxdhits->GetEntries();
101  if (nevents < minEvents) {
102  if (not forceContinueMasking) {
103  B2INFO("Not enough data: Only " << nevents << " events were collected!");
104  return c_NotEnoughData;
105  } else {
106  B2WARNING("Not enough data: Only " << nevents << " events were collected! The masking continous but the mask may be empty.");
107  }
108  }
109 
110  B2RESULT("Found total of " << nevents << " events in collected data.");
111 
112  // Get the total number of PXD hits and sensors
113  unsigned long long int nPXDHits = 0;
114  int nPXDSensors = 0;
115  for (auto sensBin = 1; sensBin <= collector_pxdhitcounts->GetXaxis()->GetNbins(); sensBin++) {
116  // The bin label is assumed to be a string representation of VxdID
117  string sensorDescr = collector_pxdhitcounts->GetXaxis()->GetBinLabel(sensBin);
118  VxdID id(sensorDescr);
119  //Increment number of sensors
120  nPXDSensors += 1;
121  // Increment number of of collected hits
122  unsigned long long int nSensorHits = collector_pxdhitcounts->GetBinContent(sensBin);
123  nPXDHits += nSensorHits;
124 
125  B2RESULT("Number of hits for sensor sensor=" << id << " is " << nSensorHits);
126  }
127 
128  // We should have enough hits in the PXD before we decide a single sensor is dead
129  unsigned long long int minPXDHits = minHits * nPXDSensors * c_nUCells * c_nVCells;
130  if (nPXDHits < minPXDHits) {
131  if (not forceContinueMasking) {
132  B2INFO("Not enough data: Only " << nPXDHits << " raw hits were collected!");
133  return c_NotEnoughData;
134  } else {
135  B2WARNING("Not enough data: Only " << nPXDHits << " raw hits were collected! The masking continous but the mask may be empty.");
136  }
137  }
138 
139  B2RESULT("Found total of " << nPXDHits << " raw hits in collected data.");
140 
141  // Check that the median number of hits is large enough
142  for (auto sensBin = 1; sensBin <= collector_pxdhitcounts->GetXaxis()->GetNbins(); sensBin++) {
143  // The bin label is assumed to be a string representation of VxdID
144  string sensorDescr = collector_pxdhitcounts->GetXaxis()->GetBinLabel(sensBin);
145  VxdID id(sensorDescr);
146 
147  m_medianOfHitsMap[id] = 0.0;
148 
149  // Get hitmap from collector
150  string name = str(format("PXD_%1%_PixelHitmap") % id.getID());
151  auto collector_pxdhitmap = getObjectPtr<TH1I>(name.c_str());
152 
153  // Check if there was data collected for this sensor
154  if (collector_pxdhitmap == nullptr) continue;
155 
156  // Compute the median number of hits to define a robust baseline for judging a channel fires too often
157  int nBins = collector_pxdhitmap->GetXaxis()->GetNbins();
158  double prob = 0.5;
159  vector<double> hitVec(nBins);
160 
161  for (auto bin = 1; bin <= nBins; bin++) {
162  hitVec[bin - 1] = (double) collector_pxdhitmap->GetBinContent(bin);
163  }
164  double medianNumberOfHits;
165  TMath::Quantiles(nBins, 1, &hitVec[0], &medianNumberOfHits, &prob, kFALSE);
166  if (medianNumberOfHits <= 0) {
167  B2WARNING("Median number of hits per senor is smaller <1. Raise median to 1 instead.");
168  medianNumberOfHits = 1;
169  } else {
170  B2RESULT("Median of hits is " << medianNumberOfHits << " for sensor " << id);
171  }
172 
173  // Keep the median of later use
174  m_medianOfHitsMap[id] = medianNumberOfHits;
175 
176  // Check median number of hits is large enough
177  if (medianNumberOfHits < minHits) {
178  if (not forceContinueMasking) {
179  B2INFO("Not enough data: Median number of his is only " << medianNumberOfHits << "!");
180  return c_NotEnoughData;
181  } else {
182  B2WARNING("Not enough data: Median number of hits is only " << medianNumberOfHits <<
183  "! The masking continous but the mask may be empty.");
184  }
185  }
186  }
187 
188  // This is the occupancy info payload for conditions DB
189  PXDOccupancyInfoPar* occupancyInfoPar = new PXDOccupancyInfoPar();
190 
191  // This is the dead pixels payload for conditions DB
192  PXDDeadPixelPar* deadPixelsPar = new PXDDeadPixelPar();
193 
194  // This is the hot pixel masking payload for conditions DB
195  PXDMaskedPixelPar* maskedPixelsPar = new PXDMaskedPixelPar();
196 
197  // Remember the number of events, helps to judge the reliability of calibrations.
198  occupancyInfoPar->setNumberOfEvents(nevents);
199 
200  // Compute the masks for all sensors
201  for (auto sensBin = 1; sensBin <= collector_pxdhitcounts->GetXaxis()->GetNbins(); sensBin++) {
202  // The bin label is assumed to be a string representation of VxdID
203  string sensorDescr = collector_pxdhitcounts->GetXaxis()->GetBinLabel(sensBin);
204  VxdID id(sensorDescr);
205 
206  // If we reach here, a sensor with no hits is deemed dead
207  int nSensorHits = collector_pxdhitcounts->GetBinContent(sensBin);
208  if (nSensorHits == 0) {
209  deadPixelsPar->maskSensor(id.getID());
210  continue;
211  }
212 
213  // Get hitmap from collector
214  string name = str(format("PXD_%1%_PixelHitmap") % id.getID());
215  auto collector_pxdhitmap = getObjectPtr<TH1I>(name.c_str());
216  if (collector_pxdhitmap == nullptr) {
217  B2WARNING("Cannot find PixelHitmap although there should be hits. This is strange!");
218  continue;
219  }
220 
221  double medianNumberOfHits = m_medianOfHitsMap[id];
222  int nBins = collector_pxdhitmap->GetXaxis()->GetNbins();
223 
224  // Dead pixel masking
225  if (medianNumberOfHits >= minHits) {
226 
227  // Bookkeeping for masking of drains and rows
228  vector<int> hitsAlongRow(c_nVCells, 0);
229  vector<int> hitsAlongDrain(c_nDrains, 0);
230 
231  // Accumulate hits along drains and rows
232  for (auto bin = 1; bin <= nBins; bin++) {
233  // Find the current pixel cell
234  int pixID = bin - 1;
235  int uCell = pixID / c_nVCells;
236  int vCell = pixID % c_nVCells;
237  int drainID = uCell * 4 + vCell % 4;
238  int nhits = collector_pxdhitmap->GetBinContent(bin);
239  hitsAlongDrain[drainID] += nhits;
240  hitsAlongRow[vCell] += nhits;
241  }
242 
243  // Dead row masking
244  B2INFO("Entering masking of dead rows ...");
245  for (auto vCell = 0; vCell < c_nVCells; vCell++) {
246  // Get number of hits per row
247  int nhits = hitsAlongRow[vCell];
248  // Mask dead row
249  if (nhits == 0) {
250  deadPixelsPar->maskRow(id.getID(), vCell);
251  B2RESULT("Dead row with vCell=" << vCell << " on sensor " << id);
252  }
253  }
254 
255  // Dead drain masking
256  B2INFO("Entering masking of dead drains ...");
257  for (auto drainID = 0; drainID < c_nDrains; drainID++) {
258  // Get number of hits per drain
259  int nhits = hitsAlongDrain[drainID];
260  // Mask dead drain
261  if (nhits == 0) {
262  deadPixelsPar->maskDrain(id.getID(), drainID);
263  B2RESULT("Dead drain line at drainID=" << drainID << " on sensor " << id);
264  }
265  }
266 
267  // Dead pixel masking
268  B2INFO("Entering masking of single dead pixels ...");
269  for (auto bin = 1; bin <= nBins; bin++) {
270  // First, we mask single pixels exceeding hit threshold
271  int nhits = collector_pxdhitmap->GetBinContent(bin);
272  int pixID = bin - 1;
273  int uCell = pixID / c_nVCells;
274  int vCell = pixID % c_nVCells;
275  int drainID = uCell * 4 + vCell % 4;
276  // Mask dead pixel
277  if (nhits == 0 && !deadPixelsPar->isDeadRow(id.getID(), vCell) && !deadPixelsPar->isDeadDrain(id.getID(), drainID)) {
278  // This pixel is dead, we have to mask it
279  deadPixelsPar->maskSinglePixel(id.getID(), pixID);
280  B2RESULT("Dead single pixel with ucell=" << uCell << ", vcell=" << vCell << " on sensor " << id);
281  }
282  }
283  }
284 
285  // Hot pixel masking
286 
287  // Bookkeeping for masking hot drains
288  vector<float> unmaskedHitsAlongDrain(c_nDrains, 0);
289  vector<int> unmaskedCellsAlongDrain(c_nDrains, 0);
290 
291  // Bookkeeping for maskign hot rows
292  vector<float> unmaskedHitsAlongRow(c_nVCells, 0);
293  vector<int> unmaskedCellsAlongRow(c_nVCells, 0);
294 
295  // Mask all single pixels exceeding medianNumberOfHits x multiplier
296  double pixelHitThr = pixelMultiplier * medianNumberOfHits;
297  B2RESULT("Pixel hit threshold is " << pixelHitThr << " for sensor " << id);
298 
299  // Mask all hot pixel for this sensor
300  for (auto bin = 1; bin <= nBins; bin++) {
301  // Find the current pixel cell
302  int pixID = bin - 1;
303  int uCell = pixID / c_nVCells;
304  int vCell = pixID % c_nVCells;
305  int drainID = uCell * 4 + vCell % 4;
306 
307  // First, we mask single pixels exceeding hit threshold
308  float nhits = collector_pxdhitmap->GetBinContent(bin);
309  bool masked = false;
310 
311  if (nhits > pixelHitThr) {
312  // This pixel is hot, we have to mask it
313  maskedPixelsPar->maskSinglePixel(id.getID(), pixID);
314  masked = true;
315  B2RESULT("Masking single pixel with ucell=" << uCell << ", vcell=" << vCell << " on sensor " << id);
316  }
317 
318  // Then we accumulate hits along u and v direction for unmasked
319  // pixels
320  if (not masked) {
321  ++unmaskedCellsAlongDrain[drainID];
322  unmaskedHitsAlongDrain[drainID] += nhits;
323  ++unmaskedCellsAlongRow[vCell];
324  unmaskedHitsAlongRow[vCell] += nhits;
325  }
326  }
327 
328  if (maskDrains) {
329  double drainHitThr = drainMultiplier * medianNumberOfHits;
330  B2RESULT("Drain hit threshold is " << drainHitThr << " for sensor " << id);
331 
332  for (auto drainID = 0; drainID < c_nDrains; drainID++) {
333  if (unmaskedCellsAlongDrain[drainID] > 0) {
334  // Compute average number of hits per drain
335  float nhits = unmaskedHitsAlongDrain[drainID] / unmaskedCellsAlongDrain[drainID];
336  // Mask residual hot drain
337  if (nhits > drainHitThr) {
338  for (auto iGate = 0; iGate < 192; iGate++) {
339  int uCell = drainID / 4;
340  int vCell = drainID % 4 + iGate * 4;
341  maskedPixelsPar->maskSinglePixel(id.getID(), uCell * c_nVCells + vCell);
342  }
343  B2RESULT("Masking drain line with drainID=" << drainID << " on sensor " << id);
344  }
345  }
346  }
347  }
348 
349  if (maskRows) {
350  double rowHitThr = rowMultiplier * medianNumberOfHits;
351  B2RESULT("Row hit threshold is " << rowHitThr << " for sensor " << id);
352 
353  for (auto vCell = 0; vCell < c_nVCells; vCell++) {
354  if (unmaskedCellsAlongRow[vCell] > 0) {
355  // Compute average number of hits per row
356  float nhits = unmaskedHitsAlongRow[vCell] / unmaskedCellsAlongRow[vCell];
357  // Mask residual hot row
358  if (nhits > rowHitThr) {
359  for (auto uCell = 0; uCell < c_nUCells; uCell++)
360  maskedPixelsPar->maskSinglePixel(id.getID(), uCell * c_nVCells + vCell);
361 
362  B2RESULT("Masking complete row with vCell=" << vCell << " on sensor " << id);
363  }
364  }
365  }
366  }
367 
368  // After the masking is done, we compute the average sensor occupancy after
369  // hot pixel masking.
370 
371  // Count all unmasked hits
372  int numberOfUnmaskedHits = 0;
373 
374  for (auto bin = 1; bin <= nBins; bin++) {
375  // Find the current pixel cell
376  int pixID = bin - 1;
377 
378  if (maskedPixelsPar->pixelOK(id.getID(), pixID)) {
379  numberOfUnmaskedHits += collector_pxdhitmap->GetBinContent(bin);
380  }
381  }
382 
383  // Compute mean occupancy before masking
384  float meanOccupancyAfterMasking = (float)numberOfUnmaskedHits / nevents / c_nVCells / c_nUCells;
385  B2RESULT("Hotpixel filtered occupancy sensor=" << id << " is " << meanOccupancyAfterMasking);
386 
387  occupancyInfoPar->setOccupancy(id.getID(), meanOccupancyAfterMasking);
388  }
389 
390  for (auto elem : maskedPixelsPar->getMaskedPixelMap()) {
391  auto id = elem.first;
392  auto singles = elem.second;
393  B2RESULT("SensorID " << VxdID(id) << " has filtered occupancy of " << occupancyInfoPar->getOccupancy(id));
394  B2RESULT("SensorID " << VxdID(id) << " has fraction of masked pixels of " << (float)singles.size() / (c_nVCells * c_nUCells));
395  }
396 
397  // Save the hot pixel mask to database. Note that this will set the database object name to the same as the collector but you
398  // are free to change it.
399  saveCalibration(maskedPixelsPar, "PXDMaskedPixelPar");
400  saveCalibration(deadPixelsPar, "PXDDeadPixelPar");
401  saveCalibration(occupancyInfoPar, "PXDOccupancyInfoPar");
402 
403  // Create debugging histo if we asked for it
405 
406  B2INFO("PXDHotPixelMask Calibration Successful");
407  return c_OK;
408 }
Belle2::PXDHotPixelMaskCalibrationAlgorithm::forceContinueMasking
bool forceContinueMasking
Force continue masking in almost empty runs instead of returning c_NotEnoughData.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:39
Belle2::PXDDeadPixelPar::maskDrain
void maskDrain(unsigned short sensorID, unsigned int drainID)
Mask single drain.
Definition: PXDDeadPixelPar.h:69
Belle2::PXDHotPixelMaskCalibrationAlgorithm::createDebugHistogram
void createDebugHistogram()
Perform debug histogram file creation.
Definition: PXDHotPixelMaskCalibrationAlgorithm.cc:39
Belle2::PXDHotPixelMaskCalibrationAlgorithm::minHits
int minHits
Minimum median number of hits per pixel needed for dead pixel masking.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:45
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::PXDHotPixelMaskCalibrationAlgorithm::drainMultiplier
float drainMultiplier
The occupancy threshold for masking hot drains is the median occupancy x drainMultiplier.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:54
Belle2::PXDHotPixelMaskCalibrationAlgorithm::m_debugHisto
int m_debugHisto
Set if a debugging histogram should be created in the algorithm output directory.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:79
Belle2::PXDHotPixelMaskCalibrationAlgorithm::m_medianOfHitsMap
std::map< VxdID, double > m_medianOfHitsMap
map of VxdID to median hits of each sensor
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:82
Belle2::PXDDeadPixelPar::isDeadRow
bool isDeadRow(unsigned short sensorID, unsigned int vCellID) const
Check whether a row is dead.
Definition: PXDDeadPixelPar.h:150
Belle2::PXDHotPixelMaskCalibrationAlgorithm::pixelMultiplier
float pixelMultiplier
The occupancy threshold for masking hot single pixels is the median occupancy x pixelMultiplier.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:48
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::getID
int getID(const std::vector< double > &breaks, double t)
get id of the time point t
Definition: calibTools.h:71
Belle2::PXDHotPixelMaskCalibrationAlgorithm::c_nDrains
const unsigned short c_nDrains
Number of drain lines of Belle II PXD sensors.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:76
Belle2::PXDHotPixelMaskCalibrationAlgorithm::calibrate
virtual EResult calibrate() override
Run algo on data.
Definition: PXDHotPixelMaskCalibrationAlgorithm.cc:93
Belle2::PXDHotPixelMaskCalibrationAlgorithm::minEvents
int minEvents
Minimum number of collected events.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:42
Belle2::PXDDeadPixelPar::maskSensor
void maskSensor(unsigned short sensorID)
Mask sensor.
Definition: PXDDeadPixelPar.h:55
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::PXDOccupancyInfoPar::setOccupancy
void setOccupancy(unsigned short sensorID, float occupancy)
Set occupancy.
Definition: PXDOccupancyInfoPar.h:57
Belle2::PXDDeadPixelPar
The payload telling which PXD pixel is dead (=Readout system does not receive signals)
Definition: PXDDeadPixelPar.h:43
Belle2::PXDHotPixelMaskCalibrationAlgorithm::maskDrains
bool maskDrains
Mask hot drain lines with too high average occupancy after single pixel masking.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::getRunListFromAllData
std::vector< Calibration::ExpRun > getRunListFromAllData() const
Get the complete list of runs from inspection of collected data.
Definition: CalibrationAlgorithm.cc:311
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::PXDOccupancyInfoPar::setNumberOfEvents
void setNumberOfEvents(int nEvents)
Set number of events used for occupancy estimation.
Definition: PXDOccupancyInfoPar.h:46
Belle2::PXDHotPixelMaskCalibrationAlgorithm::rowMultiplier
float rowMultiplier
The occupancy threshold for masking hot rows is the median occupancy x rowMultiplier.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:60
Belle2::PXDDeadPixelPar::maskSinglePixel
void maskSinglePixel(unsigned short sensorID, unsigned int pixID)
Mask single pixel.
Definition: PXDDeadPixelPar.h:115
Belle2::PXDHotPixelMaskCalibrationAlgorithm::m_file
std::shared_ptr< TFile > m_file
Pointer for TFile.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:85
Belle2::PXDMaskedPixelPar::getMaskedPixelMap
const std::unordered_map< unsigned short, MaskedSinglePixelsSet > & getMaskedPixelMap() const
Return unordered_map with all masked single pixels in PXD.
Definition: PXDMaskedPixelPar.h:87
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::CalibrationAlgorithm::getRunList
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Definition: CalibrationAlgorithm.h:276
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::PXDOccupancyInfoPar::getOccupancy
float getOccupancy(unsigned short sensorID) const
Get occupancy.
Definition: PXDOccupancyInfoPar.h:67
Belle2::PXDMaskedPixelPar::maskSinglePixel
void maskSinglePixel(unsigned short sensorID, unsigned int pixID)
Mask single pixel.
Definition: PXDMaskedPixelPar.h:49
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::PXDHotPixelMaskCalibrationAlgorithm::c_nUCells
const unsigned short c_nUCells
Number of uCells of Belle II PXD sensors.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:74
Belle2::PXDDeadPixelPar::isDeadDrain
bool isDeadDrain(unsigned short sensorID, unsigned int drainID) const
Check whether a drain is dead.
Definition: PXDDeadPixelPar.h:168
Belle2::PXDOccupancyInfoPar
The payload collecting some meta information from running the masking algorithm.
Definition: PXDOccupancyInfoPar.h:38
Belle2::PXDMaskedPixelPar::pixelOK
bool pixelOK(unsigned short sensorID, unsigned int pixID) const
Check whether a pixel on a given sensor is OK or not.
Definition: PXDMaskedPixelPar.h:72
Belle2::PXDDeadPixelPar::maskRow
void maskRow(unsigned short sensorID, unsigned int vCellID)
Mask single row.
Definition: PXDDeadPixelPar.h:92
Belle2::CalibrationAlgorithm::getPrefix
std::string getPrefix() const
Get the prefix used for getting calibration data.
Definition: CalibrationAlgorithm.h:156
Belle2::PXDHotPixelMaskCalibrationAlgorithm::c_nVCells
const unsigned short c_nVCells
Number of vCells of Belle II PXD sensors.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:72
Belle2::PXDHotPixelMaskCalibrationAlgorithm::maskRows
bool maskRows
Mask hot rows with too high average occupancy after single pixel masking.
Definition: PXDHotPixelMaskCalibrationAlgorithm.h:57
Belle2::PXDMaskedPixelPar
The payload telling which PXD pixel to mask (ignore)
Definition: PXDMaskedPixelPar.h:34