Belle II Software development
dqm.py
1#!/usr/bin/env python3
2
3
10from array import array
11from ROOT import TCanvas, TF1, TFile
12from ROOT import TGraph, TH2D, TLine
13from ROOT import kFullCircle
14import argparse
15
16nWires = [160, 160, 160, 160, 160, 160, 160, 160,
17 160, 160, 160, 160, 160, 160,
18 192, 192, 192, 192, 192, 192,
19 224, 224, 224, 224, 224, 224,
20 256, 256, 256, 256, 256, 256,
21 288, 288, 288, 288, 288, 288,
22 320, 320, 320, 320, 320, 320,
23 352, 352, 352, 352, 352, 352,
24 384, 384, 384, 384, 384, 384]
25
26nWiresSum = []
27sum = 0
28for w in nWires:
29 sum = sum + w
30 nWiresSum.append(sum)
31
32nWiresSL = [160 * 8, 160 * 6, 192 * 6, 224 * 6, 256 * 6, 288 * 6,
33 320 * 6, 352 * 6, 384 * 6]
34
35rfClock = 1.017774
36refT0 = 4905.0 * rfClock
37refMPV = 80
38
39
40class DQM():
41
42 """
43 DQM class.
44 """
45
46 def __init__(self, fname='output.root'):
47 """
48 call constructor of base class, required.
49 """
50
51 super().__init__()
52
53
54 self.f = TFile(fname)
55
56 self.type = 'adc'
57
58 self.histADC = [self.f.Get('h' + str(100000 + i))
59 for i in range(14336)]
60
61 self.histTDC = [self.f.Get('h' + str(200000 + i))
62 for i in range(14336)]
63
64 self.histHit = [self.f.Get('h' + str(400000 + i))
65 for i in range(56)]
66
67 self.lines = []
68 for h in self.histTDC:
69 y = h.GetMaximum()
70 line = TLine(refT0, 0, refT0, y)
71 line.SetLineColor(2)
72 self.lines.append(line)
73
74
75 self.lineMPVs = []
76 for h in self.histADC:
77 y = h.GetMaximum()
78 line = TLine(refMPV, 0, refMPV, y)
79 line.SetLineColor(2)
80 self.lineMPVs.append(line)
81 # self.histADCinLayer = [self.f.Get('h' + str(500000 + i))
82 # for i in range(56)]
83
84
85 self.sig = TF1("sig", "landau", 0, 200)
86 # self.sig = TF1("sig", "landau+expo(3)", 0, 200)
87 # self.sig.SetParameters(300.0,60.0,20.0,1000,-1.0)
88
89 self.ft0 = TF1("ft0", "[0]+[1]*(exp([2]*([3]-x))/(1+exp(-([4]-x)/[5])))", 2000, 4000)
90 self.ft0.SetParameters(0, 10, 0, 0, 3800, 3.)
91 self.ft0.FixParameter(0, 0.)
92
93 self.m_fitStatus = False
94
95 self.canvas = TCanvas("canvas", "canvas", 800, 800)
96
97 self.hid = 0
98
99 self.ndiv = 1
100
101 self.par = [-1.0 for i in range(14336)]
102
103 self.parADCinLayer = [-1.0 for i in range(56)]
104
105 self.parT0 = [0.0 for i in range(14336)]
106
107 self.h2 = TH2D("h2d", "MPV of all ADCs", 15000,
108 -500, 14500, 100, 0, 200)
109
110 self.h2.SetStats(0)
111
112 self.h2.GetXaxis().SetTitle('Cell ID')
113
114 self.h2.GetYaxis().SetTitle('MPV')
115
116 self.h2ADCvsLayer = TH2D("h2dl", "MPV for every Layer", 56,
117 0, 56, 100, 0, 200)
118
119 self.h2ADCvsLayer.SetStats(0)
120
121 self.h2ADCvsLayer.GetXaxis().SetTitle('Layer')
122
123 self.h2ADCvsLayer.GetYaxis().SetTitle('MPV (adc count)')
124
125 self.x = array("d", [i for i in range(14336)])
126
127 self.l_array = array("d", [i for i in range(56)])
128
129 self.graph = TGraph(len(self.x))
130
131 self.graph.SetMarkerStyle(kFullCircle)
132
133 self.graph.SetMarkerSize(0.3)
134
135
136 self.graph2 = TGraph(len(self.l_array))
137
138 self.graph2.SetMarkerStyle(kFullCircle)
139
140 self.graph2.SetMarkerSize(1.0)
141
142 self.line = [TLine(0, 0, 0, 200)]
143
144 x0 = 0
145 for i in range(9):
146 x0 = x0 + nWiresSL[i]
147 line = TLine(x0, 0, x0, 200)
148 self.line.append(line)
149 for l_i in self.line:
150 l_i.SetLineColorAlpha(8, 0.5)
151 l_i.SetLineStyle(2)
152
153 def fitADC(self, xmin=8, xmax=200):
154 """
155 Fit all ADC histgrams.
156 """
157 if self.m_fitStatus is True:
158 print("Fitting for ADC has been already applied.")
159 else:
160
161 for (i, h) in enumerate(self.histADC):
162 if h.GetEntries() > 0:
163 h.Fit(self.sig, "Q", "", xmin, xmax)
164 self.par[i] = self.sig.GetParameter(1)
165 self.m_fitStatus = True
166
167 def getHistID(self, lay=0, wire=0):
168 """
169 Get Hist ID from (layer, wire) ID.
170 Return values : hist ID (cell ID)
171 """
172 hid = nWiresSum[lay - 1] + wire if lay > 0 else wire
173 return hid
174
175 def getLayerWireID(self, hid=0):
176 """
177 Get (layer, wire) ID from cell ID.
178 Return values : [layer, wire]
179 """
180 if hid < nWiresSum[0]:
181 return [0, hid]
182 else:
183 layer = 0
184 wire = 0
185 for (i, wsum) in enumerate(nWiresSum):
186 diff = hid - wsum
187 if diff > 0 and diff < nWires[i]:
188 layer = i + 1
189 wire = diff
190 break
191 return [layer, wire]
192
193 def drawADC(self, lay=0, w=0):
194 """
195 Draw ADC histgrams w.r.t (lay,wire).
196 """
197 i = self.getHistID(lay, w)
198 for j in range(self.ndiv):
199 self.hid = i + j
200 self.canvas.cd(j + 1)
201 self.histADC[self.hid].Draw()
202 self.lineMPVs[self.hid].Draw()
203
204 self.type = 'adc'
205 self.canvas.Update()
206
207 def drawTDC(self, lay=0, w=0):
208 """
209 Draw TDC histgrams w.r.t (lay,wire).
210 """
211 i = self.getHistID(lay, w)
212 for j in range(self.ndiv):
213 self.hid = i + j
214 self.canvas.cd(j + 1)
215 self.histTDC[self.hid].Draw()
216 self.lines[self.hid].Draw()
217
218 self.canvas.Update()
219 self.type = 'tdc'
220
221 def drawHit(self, lay=0):
222 """
223 Draw Hit histgrams w.r.t layer ID.
224 """
225 for j in range(self.ndiv):
226 self.hid = lay + j
227 if self.hid > 55:
228 break
229 self.canvas.cd(j + 1)
230 self.histHit[self.hid].Draw()
231 self.canvas.Update()
232 self.type = 'hit'
233
234 def drawTDCByCell(self, i):
235 """
236 Draw TDC histgrams w.r.t cell ID (0-14336).
237 """
238 for j in range(self.ndiv):
239 self.hid = i + j
240 self.canvas.cd(j + 1)
241 self.histTDC[self.hid].Draw()
242 self.lines[self.hid].Draw()
243
244 self.type = 'tdc'
245 self.canvas.Update()
246
247 def drawADCByCell(self, i):
248 """
249 Draw ADC histgrams w.r.t cell ID (0-14336).
250 """
251 for j in range(self.ndiv):
252 self.hid = i + j
253 self.canvas.cd(j + 1)
254 self.histADC[self.hid].Draw()
255 self.lineMPVs[self.hid].Draw()
256
257 self.type = 'adc'
258 self.canvas.Update()
259
260 def next(self):
261 """
262 Show next plots.
263 """
264 self.hid = self.hid + 1
265 if self.type == 'adc':
266 self.canvas.Clear('d')
267 self.drawADCByCell(self.hid)
268 elif self.type == 'tdc':
269 self.canvas.Clear('d')
270 self.drawTDCByCell(self.hid)
271 elif self.type == 'hit':
272 self.canvas.Clear('d')
273 self.drawHit(self.hid)
274 else:
275 print('Undefined type: ' + type)
276
277 def div(self, i=1, j=1):
278 """
279 Divide Tcanvas by (i,j).
280 """
281 self.canvas.Clear()
282 self.canvas.Divide(i, j)
283 self.ndiv = i * j
284
285 def mpv(self):
286 """
287 Show the MPV of ADC distribution
288 for all wires.
289 """
290 self.canvas.Clear()
291
292 for (i, p) in enumerate(self.par):
293 self.graph.SetPoint(i, self.x[i], p)
294
295 self.h2.Draw()
296 for line in self.line:
297 line.Draw()
298 self.graph.Draw("P")
299 self.canvas.Update()
300
301 def xrange(self, xmin=-500, xmax=14500):
302 """
303 Set x range of MPV plot.
304 """
305 self.h2.GetXaxis().SetRangeUser(xmin, xmax)
306
307 def xrangeHit(self, xmin=0, xmax=350):
308 """
309 Set x range of Hit distributions.
310 """
311 for h in self.histHit:
312 h.GetXaxis().SetRangeUser(xmin, xmax)
313
314 def xrangeTDC(self, xmin=0, xmax=4000):
315 """
316 Set x range of TDC distributions.
317 """
318 for h in self.histTDC:
319 h.GetXaxis().SetRangeUser(xmin, xmax)
320
321 def xrangeADC(self, xmin=0, xmax=200):
322 """
323 Set x range of ADC distributions.
324 """
325 for h in self.histADC:
326 h.GetXaxis().SetRangeUser(xmin, xmax)
327
328 def print(self, fname='test.png'):
329 """
330 Print the current plot.
331 """
332 self.canvas.Print(fname)
333
334 def ne(self):
335 """
336 Alias to next()
337 """
338 self.next()
339
340 def dh(self, lay):
341 """
342 Alias to drawHit()
343 """
344 self.drawHit(lay)
345
346 def da(self, lay, w):
347 """
348 Alias to drawADC()
349 """
350 self.drawADC(lay, w)
351
352 def dt(self, lay, w):
353 """
354 Alias to drawTDC()
355 """
356 self.drawTDC(lay, w)
357
358
359if __name__ == "__main__":
360 # Make the parameters accessible form the outside.
361
362 parser = argparse.ArgumentParser()
363 parser.add_argument('input', help='Input file to be processed (unpacked CDC data).')
364 args = parser.parse_args()
365 dqm = DQM(args.input)
Definition: dqm.py:40
def print(self, fname='test.png')
Definition: dqm.py:328
def div(self, i=1, j=1)
Definition: dqm.py:277
lines
Line for guide of TDC T0.
Definition: dqm.py:67
def xrangeADC(self, xmin=0, xmax=200)
Definition: dqm.py:321
def dh(self, lay)
Definition: dqm.py:340
def xrange(self, xmin=-500, xmax=14500)
Definition: dqm.py:301
def drawHit(self, lay=0)
Definition: dqm.py:221
line
Marker style.
Definition: dqm.py:142
def fitADC(self, xmin=8, xmax=200)
Definition: dqm.py:153
def da(self, lay, w)
Definition: dqm.py:346
def xrangeTDC(self, xmin=0, xmax=4000)
Definition: dqm.py:314
f
input file name
Definition: dqm.py:54
parT0
T0 values.
Definition: dqm.py:105
histHit
Hit distribution histograms.
Definition: dqm.py:64
histADC
ADC histograms.
Definition: dqm.py:58
m_fitStatus
Fit status True if fitting has been finished.
Definition: dqm.py:93
type
plot type (adc/tdc/hit)
Definition: dqm.py:56
h2ADCvsLayer
Set stats.
Definition: dqm.py:116
histTDC
TDC histograms.
Definition: dqm.py:61
par
MPV values.
Definition: dqm.py:101
graph
TGraph.
Definition: dqm.py:129
x
Set stats.
Definition: dqm.py:125
def getHistID(self, lay=0, wire=0)
Definition: dqm.py:167
sig
Fitting function for ADC histogram.
Definition: dqm.py:85
def next(self)
Definition: dqm.py:260
def ne(self)
Definition: dqm.py:334
hid
histogram ID
Definition: dqm.py:97
graph2
Marker style.
Definition: dqm.py:136
ft0
Fitting function for TDC histogram.
Definition: dqm.py:89
def mpv(self)
Definition: dqm.py:285
def drawTDCByCell(self, i)
Definition: dqm.py:234
def xrangeHit(self, xmin=0, xmax=350)
Definition: dqm.py:307
lineMPVs
Line for guide of ADC MPV.
Definition: dqm.py:75
def drawTDC(self, lay=0, w=0)
Definition: dqm.py:207
h2
histogram
Definition: dqm.py:107
def __init__(self, fname='output.root')
Definition: dqm.py:46
parADCinLayer
MPV values w.r.t.
Definition: dqm.py:103
ndiv
Number of division for canvas.
Definition: dqm.py:99
def drawADCByCell(self, i)
Definition: dqm.py:247
def drawADC(self, lay=0, w=0)
Definition: dqm.py:193
canvas
canvas
Definition: dqm.py:95
def getLayerWireID(self, hid=0)
Definition: dqm.py:175
l_array
array layer
Definition: dqm.py:127