Belle II Software  release-06-02-00
plotHotPixelMasks.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # Creates overview plots for hotpixel calibrations
13 #
14 # At first, you can extract the hotpixel calibration payloads from a localdb/centraldb using the tool
15 #
16 # b2conditionsdb-extract --exp 3 --runs 0-5614 --tag Calibration_Offline_Development --output hot_payloads.root PXDMaskedPixelPar
17 #
18 # Secondly, execute the script as
19 #
20 # basf2 plotHotPixelMasks.py
21 #
22 # basf2 plotHotPixelMasks.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 hotpixel 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('hotpixel_histos.root', 'RECREATE')
50 histofile.cd()
51 histofile.mkdir("maps")
52 
53 # Open file with extracted payloads
54 rfile = ROOT.TFile("hot_payloads.root", "READ")
55 b2.conditions = rfile.Get("conditions")
56 
57 hotpixel_table = dict()
58 for sensorID in sensor_list:
59  hotpixel_table[sensorID.getID()] = list()
60 run_list = list()
61 
62 for condition in b2.conditions:
63  if condition.PXDMaskedPixelPar_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  hotpixelmap = condition.PXDMaskedPixelPar.getMaskedPixelMap()
73 
74  layer = sensorID.getLayerNumber()
75  ladder = sensorID.getLadderNumber()
76  sensor = sensorID.getSensorNumber()
77 
78  counter = len(hotpixelmap[sensorID.getID()])
79 
80  if args.maps:
81  name = "MaskedPixels_{:d}_{:d}_{:d}_run_{:d}".format(layer, ladder, sensor, condition.run)
82  title = "Masked Pixels Sensor={:d}.{:d}.{:d} run={:d}".format(layer, ladder, sensor, condition.run)
83  hot_map = ROOT.TH2F(name, title, nUCells, 0, nUCells, nVCells, 0, nVCells)
84  hot_map.GetXaxis().SetTitle("uCell")
85  hot_map.GetYaxis().SetTitle("vCell")
86  hot_map.GetZaxis().SetTitle("isMasked")
87  hot_map.SetStats(0)
88 
89  for uCell in range(nUCells):
90  for vCell in range(nVCells):
91  pixID = uCell * nVCells + vCell
92  isMasked = not condition.PXDMaskedPixelPar.pixelOK(sensorID.getID(), pixID)
93  hot_map.SetBinContent(int(uCell + 1), int(vCell + 1), isMasked)
94 
95  histofile.cd("maps")
96  hot_map.Write()
97 
98  hotfraction = counter / (nUCells * nVCells)
99  hotpixel_table[sensorID.getID()].append(hotfraction)
100 
101 
102 histofile.cd()
103 c = ROOT.TCanvas('hotpixels_vs_runno', 'Hotpixel evolution vs. run number', 200, 10, 700, 500)
104 c.SetGrid()
105 
106 for sensorID in sensor_list:
107 
108  n = len(run_list)
109  x, y = array('d'), array('d')
110  for value in hotpixel_table[sensorID.getID()]:
111  y.append(value)
112  for run in run_list:
113  x.append(run)
114 
115  gr = ROOT.TGraph(n, x, y)
116  gr.SetLineColor(ROOT.kBlue)
117  gr.SetLineWidth(4)
118  gr.SetName("graph_{}".format(sensorID.getID()))
119  gr.SetMarkerColor(ROOT.kBlue)
120  gr.SetMarkerStyle(21)
121  gr.SetTitle('Hotpixel evolution Sensor={}'.format(sensorID))
122  gr.GetXaxis().SetTitle('run number')
123  gr.GetYaxis().SetTitle('hotpixel fraction')
124  gr.GetYaxis().SetRangeUser(0.0, 1.0)
125  gr.Draw('AP')
126 
127  c.Update()
128  c.Modified()
129  c.Print('hotpixel_vs_runno_{}.png'.format(sensorID.getID()))
130  c.Write()
131 
132 
133 rfile.Close()
134 histofile.Close()
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33