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