Belle II Software  release-06-00-14
plotDeadPixelMasks.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # Creates overview plots for deadpixel calibrations
13 #
14 # At first, you can extract the deadpixel calibration payloads from a localdb/centraldb using the tool
15 #
16 # b2conditionsdb-extract --exp 3 --runs 0-5614 --tag Calibration_Offline_Development --output dead_payloads.root PXDDeadPixelPar
17 #
18 # Secondly, execute the script as
19 #
20 # basf2 plotDeadPixelMasks.py
21 #
22 # basf2 plotDeadPixelMasks.py -- --maps
23 
24 
25 import basf2 as b2
26 import ROOT
27 from ROOT import Belle2
28 from array import array
29 
30 import argparse
31 parser = argparse.ArgumentParser(description="Plot dead pixel maps")
32 parser.add_argument('--maps', dest='maps', action="store_true", help='Create maps from payloads. This can be slow!')
33 args = parser.parse_args()
34 
35 
36 sensor_list = [
37  Belle2.VxdID("1.1.1"), Belle2.VxdID("1.1.2"),
38  Belle2.VxdID("1.2.1"), Belle2.VxdID("1.2.2"),
39  Belle2.VxdID("1.3.1"), Belle2.VxdID("1.3.2"),
40  Belle2.VxdID("1.4.1"), Belle2.VxdID("1.4.2"),
41  Belle2.VxdID("1.5.1"), Belle2.VxdID("1.5.2"),
42  Belle2.VxdID("1.6.1"), Belle2.VxdID("1.6.2"),
43  Belle2.VxdID("1.7.1"), Belle2.VxdID("1.7.2"),
44  Belle2.VxdID("1.8.1"), Belle2.VxdID("1.8.2"),
45  Belle2.VxdID("2.4.1"), Belle2.VxdID("2.4.2"),
46  Belle2.VxdID("2.5.1"), Belle2.VxdID("2.5.2")]
47 
48 # Create output file wíth histos and plots
49 histofile = ROOT.TFile('deadpixel_histos.root', 'RECREATE')
50 histofile.cd()
51 histofile.mkdir("maps")
52 
53 # Open file with extracted payloads
54 rfile = ROOT.TFile("dead_payloads.root", "READ")
55 b2.conditions = rfile.Get("conditions")
56 
57 deadpixel_table = dict()
58 for sensorID in sensor_list:
59  deadpixel_table[sensorID.getID()] = list()
60 run_list = list()
61 
62 for condition in b2.conditions:
63  if condition.PXDDeadPixelPar_valid:
64  print("Starting on run {}".format(condition.run))
65  run_list.append(condition.run)
66 
67  for sensorID in sensor_list:
68 
69  nUCells = 250
70  nVCells = 768
71 
72  deadsensormap = condition.PXDDeadPixelPar.getDeadSensorMap()
73  deaddrainmap = condition.PXDDeadPixelPar.getDeadDrainMap()
74  deadrowmap = condition.PXDDeadPixelPar.getDeadRowMap()
75  deadsinglesmap = condition.PXDDeadPixelPar.getDeadSinglePixelMap()
76 
77  layer = sensorID.getLayerNumber()
78  ladder = sensorID.getLadderNumber()
79  sensor = sensorID.getSensorNumber()
80 
81  # Every dead drain counts for 192 dead pixels
82  counter = len(deaddrainmap[sensorID.getID()]) * 192
83 
84  # Every dead row counts for 250 dead pixels
85  # This can lead to double counting of dead pixel from dead drains
86  # The double counting can be avoided using --maps option.
87  counter += len(deadrowmap[sensorID.getID()]) * 250
88 
89  # Every dead row counts for 250 dead pixels
90  counter += len(deadsinglesmap[sensorID.getID()])
91 
92  if args.maps:
93  counter = 0
94 
95  name = "DeadPixels_{:d}_{:d}_{:d}_run_{:d}".format(layer, ladder, sensor, condition.run)
96  title = "Dead Pixels Sensor={:d}.{:d}.{:d} run={:d}".format(layer, ladder, sensor, condition.run)
97  dead_map = ROOT.TH2F(name, title, nUCells, 0, nUCells, nVCells, 0, nVCells)
98  dead_map.GetXaxis().SetTitle("uCell")
99  dead_map.GetYaxis().SetTitle("vCell")
100  dead_map.GetZaxis().SetTitle("isDead")
101  dead_map.SetStats(0)
102 
103  for uCell in range(nUCells):
104  for vCell in range(nVCells):
105  pixID = uCell * nVCells + vCell
106  isDeadSinglePixel = condition.PXDDeadPixelPar.isDeadSinglePixel(sensorID.getID(), pixID)
107  isDeadRow = condition.PXDDeadPixelPar.isDeadRow(sensorID.getID(), vCell)
108  isDeadDrain = condition.PXDDeadPixelPar.isDeadDrain(sensorID.getID(), uCell * 4 + vCell % 4)
109  isDead = isDeadSinglePixel or isDeadRow or isDeadDrain
110 
111  if isDead:
112  counter += 1.0
113  dead_map.SetBinContent(int(uCell + 1), int(vCell + 1), isDead)
114 
115  histofile.cd("maps")
116  dead_map.Write()
117 
118  deadfraction = counter / (nUCells * nVCells)
119  deadpixel_table[sensorID.getID()].append(deadfraction)
120 
121 
122 histofile.cd()
123 c = ROOT.TCanvas('dead_vs_runno', 'Deadpixel evolution vs. run number', 200, 10, 700, 500)
124 c.SetGrid()
125 
126 for sensorID in sensor_list:
127 
128  n = len(run_list)
129  x, y = array('d'), array('d')
130  for value in deadpixel_table[sensorID.getID()]:
131  y.append(value)
132  for run in run_list:
133  x.append(run)
134 
135  gr = ROOT.TGraph(n, x, y)
136  gr.SetLineColor(ROOT.kBlue)
137  gr.SetLineWidth(4)
138  gr.SetName("graph_{}".format(sensorID.getID()))
139  gr.SetMarkerColor(ROOT.kBlue)
140  gr.SetMarkerStyle(21)
141  gr.SetTitle('Deadpixel evolution Sensor={}'.format(sensorID))
142  gr.GetXaxis().SetTitle('run number')
143  gr.GetYaxis().SetTitle('dead fraction')
144  gr.GetYaxis().SetRangeUser(0.0, 1.0)
145  gr.Draw('AP')
146 
147  c.Update()
148  c.Modified()
149  c.Print('deadpixels_vs_runno_{}.png'.format(sensorID.getID()))
150  c.Write()
151 
152 
153 rfile.Close()
154 histofile.Close()
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33