Belle II Software  release-08-01-10
monitor_module.py
1 # !/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 import os
13 import basf2 as b2
14 from ROOT import Belle2, TCanvas, TColor, TH1F, TH2F, gStyle
15 import numpy
16 from math import pi
17 
18 
19 # make the monitoring plots
20 file_type = 'png'
21 
22 color1501 = TColor(1501, 35 / 255., 55 / 255., 59 / 255.)
23 color1502 = TColor(1502, 99 / 255., 120 / 255., 173 / 255.)
24 color1503 = TColor(1503, 235 / 255., 129 / 255., 27 / 255.)
25 color1504 = TColor(1504, 99 / 255., 173 / 255., 132 / 255.)
26 color1505 = TColor(1505, 123 / 255., 111 / 255., 37 / 255.)
27 
28 nWiresInLayer = [160, 192, 256, 320, 384]
29 
30 
31 def set_style(hist, color=1503, marker=20):
32  hist.SetLineColor(color)
33  hist.SetLineWidth(2)
34  hist.SetFillStyle(0)
35  hist.SetMarkerStyle(marker)
36  hist.SetMarkerSize(2)
37  hist.SetMarkerColor(color)
38 
39 
40 hall = TH2F('all', 'all', 5, 0, 5, 384, 0, 384)
41 hhit = TH2F('hit', 'hit', 5, 0, 5, 384, 0, 384)
42 hno_ghost = TH2F('ng', 'hit', 5, 0, 5, 384, 0, 384)
43 hghost = TH2F('ghost', 'hit', 5, 0, 5, 384, 0, 384)
44 
45 hlr = TH1F('lr', 'left right', 4, 0, 4)
46 hpri = TH1F('pri', 'priority position', 4, 0, 4)
47 hpritime = TH1F('pritime', 'priority time', 512, 0, 512)
48 hftime = TH1F('found time', 'found time', 48, -35, 13)
49 
50 used = numpy.zeros(2336, numpy.int)
51 in_tsim = numpy.zeros(2336, numpy.int)
52 
53 
54 class Monitor(b2.Module):
55  """
56  Module to make some monitor plots for TSF
57  """
58 
59  def initialize(self):
60  """
61  Initilization
62  """
63 
64  self.event_infoevent_info = Belle2.PyStoreObj('EventMetaData')
65 
66  self.tshitstshits = Belle2.PyStoreArray('CDCTriggerSegmentHits')
67 
68  self.simhitssimhits = Belle2.PyStoreArray('TSimSegmentHits')
69 
70  self.first_runfirst_run = self.event_infoevent_info.getRun()
71 
72  self.neventsnevents = 0
73 
74  def event(self):
75  """
76  event function
77  """
78  self.neventsnevents += 1
79  used.fill(0)
80  in_tsim.fill(0)
81  for simhit in self.simhitssimhits:
82  if simhit.getISuperLayer() % 2 == 0:
83  iax, iw = simhit.getISuperLayer() // 2, simhit.getIWire()
84  if simhit.getPriorityPosition() == 1:
85  iw += 1
86  if iw == nWiresInLayer[iax]:
87  iw = 0
88  hall.Fill(iax, iw)
89  in_tsim[simhit.getSegmentID()] = 1
90  for hit in self.tshitstshits:
91  ind = hit.getSegmentID()
92  # only fill the histogram if the same ID is not found in the event
93  # In other words, it can't be the same hit from another 2D
94  # even if the unpacker does not eliminate repeated hits.
95  iax, iw = hit.getISuperLayer() // 2, hit.getIWire()
96  if hit.getPriorityPosition() == 1:
97  iw += 1
98  if iw == nWiresInLayer[iax]:
99  iw = 0
100  if not used[ind]:
101  used[ind] = 1
102  hhit.Fill(iax, iw)
103  if in_tsim[ind]:
104  hno_ghost.Fill(iax, iw)
105  if not in_tsim[ind]:
106  hghost.Fill(iax, iw)
107  hlr.Fill(hit.getLeftRight())
108  hpri.Fill(hit.getPriorityPosition())
109  hpritime.Fill(hit.priorityTime())
110  hftime.Fill(hit.foundTime())
111 
112  def terminate(self):
113  """
114  termination
115  """
116  if self.neventsnevents == 0:
117  b2.B2WARNING("The monitor module is never called.\n" +
118  "There seems to be no CDC Trigger data at all!")
119  return
120  elif hhit.GetEntries() == 0:
121  b2.B2WARNING("No recorded TS hits at all!")
122  return
123  elif hall.GetEntries() == 0:
124  b2.B2ERROR("No simulated TS hits at all!")
125  return
126  if not os.path.exists('monitor_plots'):
127  os.mkdir('monitor_plots')
128  can = TCanvas('can2', 'can2', 800, 700)
129  gStyle.SetOptStat(0)
130 
131  nts = [160, 192, 256, 320, 384]
132 
133  slall = [[]] * 5
134  slhit = [[]] * 5
135  slghost = [[]] * 5
136  quos = [[]] * 5
137  slng = [[]] * 5
138  quong = [[]] * 5
139 
140  # 1 bin corresponds to 8 TS (2 bins for a merger unit)
141  width = 8
142  for iax in range(5):
143  ibin = iax + 1
144  p = hall.ProjectionY('all{}'.format(iax), ibin, ibin)
145  p.SetBins(nts[iax], 0, 2 * pi)
146  p.Rebin(width)
147  slall[iax] = p
148 
149  q = hhit.ProjectionY('hit{}'.format(iax), ibin, ibin)
150  q.SetBins(nts[iax], 0, 2 * pi)
151  q.Rebin(width)
152  slhit[iax] = q
153 
154  g = hghost.ProjectionY('ghost{}'.format(iax), ibin, ibin)
155  g.SetBins(nts[iax], 0, 2 * pi)
156  g.Rebin(width)
157  slghost[iax] = g
158 
159  s = hno_ghost.ProjectionY('ng{}'.format(iax), ibin, ibin)
160  s.SetBins(nts[iax], 0, 2 * pi)
161  s.Rebin(width)
162  slng[iax] = s
163 
164  r = g.Clone('ghost{}'.format(iax))
165  r2 = s.Clone('quong{}'.format(iax))
166  g.Sumw2()
167  p.Sumw2()
168  s.Sumw2()
169  # get the ratio of (found hit / simulated hit) with binomial error
170  r.Divide(g, p, 1.0, 1.0, 'B')
171  quos[iax] = r
172  r2.Divide(s, p, 1.0, 1.0, 'B')
173  quong[iax] = r2
174 
175  for hits, name in [(slhit, 'data TS'),
176  (slall, 'TSIM')]:
177  for i in range(5):
178  h = hits[i]
179  h.SetTitle(f'SL{i * 2}')
180  set_style(h, 1501 + i, 20 + i)
181  if name in ['data TS', 'TSIM']:
182  try:
183  h.Scale(1 / h.Integral('width'))
184  except ZeroDivisionError:
185  b2.B2WARNING(f'Not a single hit in SL{i * 2}!')
186  continue
187  height = max([g.GetMaximum() for g in hits])
188  for h in hits:
189  h.SetMaximum(1.1 * height)
190  h.SetMinimum(0)
191  options = 'L M E0'
192  if h.GetTitle() != 'SL0':
193  options += ' same'
194  h.Draw(options)
195  can.BuildLegend(.85, .76, .95, .95)
196  hits[0].SetTitle(name + ' hit distribution in run {}; #phi (rad)'.format(
197  self.first_runfirst_run))
198  can.SaveAs('monitor_plots/' + name.split()[0] +
199  '_ts_hits_{:05d}.{}'.format(self.first_runfirst_run, file_type))
200 
201  for ratio, name in [(quos, 'ghost rate'), (quong, 'efficiency (w.r.t. fast TSIM)')]:
202  upp = max([g.GetMaximum() for g in ratio])
203  low = min([g.GetMinimum() for g in ratio])
204  for i in range(5):
205  h = ratio[i]
206  h.SetTitle('SL{}'.format(2 * i))
207  h.SetMaximum(1.1 * upp)
208  h.SetMinimum(0.9 * low)
209  set_style(h, 1501 + i, 20 + i)
210  options = 'e0'
211  if i != 0:
212  options += ' same'
213  h.Draw(options)
214  can.BuildLegend(.91, .76, .98, .95)
215  ratio[0].SetTitle('TSF ' + name + ' in run {};#phi (rad)'.format(
216  self.first_runfirst_run))
217  file_name = name.split()[0]
218  can.SaveAs('monitor_plots/ts_{}_{:05d}.{}'.format(file_name, self.first_runfirst_run, file_type))
219 
220  # ghost distribution
221  for h in [hlr, hpri, hpritime, hftime]:
222  set_style(h, 1503)
223  h.Draw()
224  name = h.GetTitle().replace(' ', '_')
225  can.SaveAs('monitor_plots/ghost_{}_{:05d}.{}'.format(name, self.first_runfirst_run, file_type))
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
nevents
Number of events.
static ExpRun getRun(std::map< ExpRun, std::pair< double, double >> runs, double t)
Get exp number + run number from time.
Definition: Splitter.cc:262