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