Belle II Software  release-06-02-00
EventInspectorLookback.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # Purpose:
13 # basf module to histogram useful values in RawKLM and BKLMHit2d data-objects in an SROOT file.
14 #
15 
16 import basf2
17 import bklmDB
18 import math
19 import array
20 import ROOT
21 from ROOT import Belle2
22 
23 
24 class EventInspectorLookback(basf2.Module):
25  """Analyze RPC lookback-window parameter settings, fill histograms"""
26 
27 
28  BKLM_ID = 0x07000000
29 
30  EKLM_ID = 0x08000000
31 
32  BKLM_STRIP_BIT = 0
33 
34  BKLM_PLANE_BIT = 6
35 
36  BKLM_LAYER_BIT = 7
37 
38  BKLM_SECTOR_BIT = 11
39 
40  BKLM_SECTION_BIT = 14
41 
42  BKLM_MAXSTRIP_BIT = 15
43 
44  BKLM_STRIP_MASK = 0x3f
45 
46  BKLM_PLANE_MASK = (1 << BKLM_PLANE_BIT)
47 
48  BKLM_LAYER_MASK = (15 << BKLM_LAYER_BIT)
49 
50  BKLM_SECTOR_MASK = (7 << BKLM_SECTOR_BIT)
51 
52  BKLM_SECTION_MASK = (1 << BKLM_SECTION_BIT)
53 
54  BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
55 
56  BKLM_MODULEID_MASK = (BKLM_SECTION_MASK | BKLM_SECTOR_MASK | BKLM_LAYER_MASK)
57 
58  def __init__(self, exp, run, histName, pdfName, mode, window):
59  """Constructor
60 
61  Arguments:
62  exp (str): formatted experiment number
63  run (str): formatter run number
64  histName (str): path name of the output histogram ROOT file
65  pdfName (str): path name of the output histogram PDF file
66  mode (int): specifies the lookback-window mode
67  0: coarse window start values
68  1: coarse window width values
69  2: fine window start values
70  3: fine window width values
71  window (int, int, int): specifies the lookback-window min, max and step values
72  """
73  super().__init__()
74 
75  self.expexp = exp
76 
77  self.runrun = run
78 
79  self.histNamehistName = histName
80 
81  self.pdfNamepdfName = pdfName
82 
83  windowModes = {0: "coarse start", 1: "coarse width", 2: "fine start", 3: "fine width"}
84 
85  self.windowModewindowMode = windowModes[mode]
86 
87  self.windowMinValuewindowMinValue = window[0]
88 
89  self.windowMaxValuewindowMaxValue = window[1]
90 
91  self.windowStepValuewindowStepValue = window[2]
92  print("Mode = {0} start = {1} end = {2} step = {3}".format(mode, window[0], window[1], window[2]))
93 
94  def initialize(self):
95  """Handle job initialization: fill the mapping database, create histograms"""
96 
97  expRun = 'e{0:02d}r{1}: '.format(int(self.expexp), int(self.runrun))
98 
99 
100  self.electIdToModuleIdelectIdToModuleId = bklmDB.fillDB()
101 
102  self.sectorFBToDCsectorFBToDC = [11, 15, 2, 6, 10, 14, 3, 7, 9, 13, 0, 4, 8, 12, 1, 5]
103 
104  self.dcToSectorFBdcToSectorFB = [10, 14, 2, 6, 11, 15, 3, 7, 12, 8, 4, 0, 13, 9, 5, 1]
105 
106  self.histogramFilehistogramFile = ROOT.TFile.Open(self.histNamehistName, "RECREATE")
107  # All histograms/scatterplots in the output file will show '# of events' only
108  ROOT.gStyle.SetOptStat(10)
109 
110  # create the rawKLM histograms
111 
112 
113  self.hist_mappedRPCTimeCalhist_mappedRPCTimeCal = ROOT.TH1F(
114  'mappedRPCTimeCal', expRun + 'RPC time distribution;t - t(trigger) (ns)', 256, -0.5, 1023.5)
115 
116  self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow = {}
117 
118  self.dict_nRawKLMsdict_nRawKLMs = {}
119  for window in range(self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+1, self.windowStepValuewindowStepValue):
120  label = "mappedRPCTimeCalByWindow{0:04x}".format(window)
121  title = "{0}RPC time distribution for lookback-window {1} = {2:04x}".format(
122  expRun, self.windowModewindowMode, window) + ";t - t(trigger) (ns)"
123  self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow[window] = ROOT.TH1F(label, title, 256, -0.5, 1023.5)
124  self.dict_nRawKLMsdict_nRawKLMs[window] = 0
125 
126  self.hist_mappedRPCTimeCalByWindowhist_mappedRPCTimeCalByWindow = self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow[self.windowMinValuewindowMinValue]
127 
128  # Create the BKLMHit2d-related histograms
129 
130 
131  self.hist_occupancyForwardXYhist_occupancyForwardXY = ROOT.TH2F('occupancyForwardXY',
132  expRun + 'Forward xy occupancy;x(cm);y(cm)',
133  230, -345.0, 345.0, 230, -345.0, 345.0)
134 
135  self.hist_occupancyBackwardXYhist_occupancyBackwardXY = ROOT.TH2F('occupancyBackwardXY',
136  expRun + 'Backward xy occupancy;x(cm);y(cm)',
137  230, -345.0, 345.0, 230, -345.0, 345.0)
138 
139  self.hist_occupancyXYByWindowhist_occupancyXYByWindow = self.dict_occupancyXYByWindowdict_occupancyXYByWindow[self.windowMinValuewindowMinValue]
140 
141 
142  self.dict_occupancyXYByWindowdict_occupancyXYByWindow = {}
143 
144  self.dict_nHit2dsdict_nHit2ds = {}
145  for window in range(self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+1, self.windowStepValuewindowStepValue):
146  label = "occupancyXYByWindow{0:04x}".format(window)
147  title = "{0}Forward xy occupancy for lookback-window {1} = {2:04x};x(cm);y(cm)".format(expRun, self.windowModewindowMode, window)
148  self.dict_occupancyXYByWindowdict_occupancyXYByWindow[window] = ROOT.TH2F(label, title, 230, -345.0, 345.0, 230, -345.0, 345.0)
149  self.dict_nHit2dsdict_nHit2ds[window] = 0
150 
151 
152  self.dict_nEventsdict_nEvents = {}
153  for window in range(self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+1, self.windowStepValuewindowStepValue):
154  self.dict_nEventsdict_nEvents[window] = 0
155 
156 
157  self.windowValuewindowValue = -1
158 
159  def terminate(self):
160  """Handle job termination: draw histograms, close output files"""
161 
162  canvas = ROOT.TCanvas("canvas", self.pdfNamepdfName, 1600, 1600)
163  title = '{0}['.format(self.pdfNamepdfName)
164  canvas.SaveAs(title)
165  canvas.Clear()
166  canvas.Divide(1, 1)
167  canvas.GetPad(0).SetGrid(1, 1)
168  canvas.GetPad(0).Update()
169  nWindows = int((self.windowMaxValuewindowMaxValue - self.windowMinValuewindowMinValue) / self.windowStepValuewindowStepValue) + 1
170  hist_nEvents = ROOT.TH1F('nEvents',
171  'Number of events;Lookback-window {0} value'.format(self.windowModewindowMode),
172  nWindows, self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+self.windowStepValuewindowStepValue)
173  hist_nEvents.SetStats(False)
174  hist_nEvents.SetMinimum(0)
175  values = array.array('d', [0]) # dummy vaue for the histogram's underflow bin
176  for key in self.dict_nEventsdict_nEvents:
177  values.append(self.dict_nEventsdict_nEvents[key])
178  values.append(0) # dummy value for the histogram's overflow bin
179  hist_nEvents.SetContent(values)
180  hist_nEvents.Draw("HIST")
181  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(hist_nEvents.GetName()))
182  hist_nRawKLMs = ROOT.TH1F('nRawKLMs',
183  'Mean number of RawKLM hits per event;Lookback-window {0} value'.format(self.windowModewindowMode),
184  nWindows, self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+self.windowStepValuewindowStepValue)
185  hist_nRawKLMs.SetStats(False)
186  hist_nRawKLMs.SetMinimum(0)
187  ratios = array.array('d', [0]) # dummy ratio for the histogram's underflow bin
188  errors = array.array('d', [0]) # dummy error for the histogram's underflow bin
189  for key in self.dict_nRawKLMsdict_nRawKLMs:
190  numerator = self.dict_nRawKLMsdict_nRawKLMs[key]
191  denominator = float(self.dict_nEventsdict_nEvents[key])
192  if denominator > 0:
193  ratio = numerator / denominator
194  ratios.append(ratio)
195  errors.append(math.sqrt(ratio * (ratio + 1.0) / denominator)) # avoid 1/numerator
196  else:
197  ratios.append(0)
198  errors.append(0)
199  ratios.append(0) # dummy ratio for the histogram's overflow bin
200  errors.append(0) # dummy error for the histogram's overflow bin
201  hist_nRawKLMs.SetContent(ratios)
202  hist_nRawKLMs.SetError(errors)
203  hist_nRawKLMs.Draw("E0 X0 L")
204  hist_nRawKLMs.Draw("HIST SAME")
205  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(hist_nRawKLMs.GetName()))
206  hist_nHit2ds = ROOT.TH1F('nBKLMHit2ds',
207  'Mean number of BKLMHit2ds per event;Lookback-window {0} value'.format(self.windowModewindowMode),
208  nWindows, self.windowMinValuewindowMinValue, self.windowMaxValuewindowMaxValue+self.windowStepValuewindowStepValue)
209  hist_nHit2ds.SetStats(False)
210  hist_nHit2ds.SetMinimum(0)
211  ratios = array.array('d', [0]) # dummy ratio for the histogram's underflow bin
212  errors = array.array('d', [0]) # dummy error for the histogram's underflow bin
213  for key in self.dict_nHit2dsdict_nHit2ds:
214  numerator = self.dict_nHit2dsdict_nHit2ds[key]
215  denominator = float(self.dict_nEventsdict_nEvents[key])
216  if denominator > 0:
217  ratio = numerator / denominator
218  ratios.append(ratio)
219  errors.append(math.sqrt(ratio * (ratio + 1.0) / denominator)) # avoid 1/numerator
220  else:
221  ratios.append(0)
222  errors.append(0)
223  ratios.append(0) # dummy ratio for the histogram's overflow bin
224  errors.append(0) # dummy error for the histogram's overflow bin
225  hist_nHit2ds.SetContent(ratios)
226  hist_nHit2ds.SetError(errors)
227  hist_nHit2ds.Draw("E0 X0 L")
228  hist_nHit2ds.Draw("HIST SAME")
229  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(hist_nHit2ds.GetName()))
230 
231  self.hist_mappedRPCTimeCalhist_mappedRPCTimeCal.Draw()
232  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(self.hist_mappedRPCTimeCalhist_mappedRPCTimeCal.GetName()))
233  for key in self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow:
234  theHist = self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow[key]
235  theHist.Draw()
236  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(theHist.GetName()))
237 
238  self.hist_occupancyBackwardXYhist_occupancyBackwardXY.Draw("colz")
239  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(self.hist_occupancyBackwardXYhist_occupancyBackwardXY.GetName()))
240  self.hist_occupancyForwardXYhist_occupancyForwardXY.Draw("colz")
241  canvas.Print(self.pdfNamepdfName, "Title:{0}".format(self.hist_occupancyForwardXYhist_occupancyForwardXY.GetName()))
242  for key in self.dict_occupancyXYByWindowdict_occupancyXYByWindow:
243  theHist = self.dict_occupancyXYByWindowdict_occupancyXYByWindow[key]
244  theHist.Draw("colz")
245  lastTitle = "Title:{0}".format(theHist.GetName())
246  canvas.Print(self.pdfNamepdfName, lastTitle)
247  pdfNameLast = '{0}]'.format(self.pdfNamepdfName)
248  canvas.Print(pdfNameLast, lastTitle)
249  self.histogramFilehistogramFile.Write()
250  self.histogramFilehistogramFile.Close()
251  print('Goodbye')
252 
253  def beginRun(self):
254  """Handle begin of run: print diagnostic message"""
255  EventMetaData = Belle2.PyStoreObj('EventMetaData')
256  print('beginRun', EventMetaData.getRun())
257 
258  def endRun(self):
259  """Handle end of run: print diagnostic message"""
260  EventMetaData = Belle2.PyStoreObj('EventMetaData')
261  print('endRun', EventMetaData.getRun())
262 
263  def event(self):
264  """Process one event: fill histograms"""
265 
266  EventMetaData = Belle2.PyStoreObj('EventMetaData')
267  event = EventMetaData.getEvent()
268  rawklms = Belle2.PyStoreArray('RawKLMs')
269  hit2ds = Belle2.PyStoreArray('BKLMHit2ds')
270 
271  # Process the RawKLMs
272 
273  for copper in range(0, len(rawklms)):
274  rawklm = rawklms[copper]
275  if rawklm.GetNumEntries() != 1:
276  print('##0 Event', event, 'copper', copper, ' getNumEntries=', rawklm.GetNumEntries())
277  continue
278  nodeID = rawklm.GetNodeID(0) - self.BKLM_IDBKLM_ID
279  if nodeID >= self.EKLM_IDEKLM_ID - self.BKLM_IDBKLM_ID:
280  nodeID = nodeID - (self.EKLM_IDEKLM_ID - self.BKLM_IDBKLM_ID) + 4
281  if (nodeID < 0) or (nodeID > 4): # examine BKLM nodes only
282  continue
283  trigCtime = (rawklm.GetTTCtime(0) & 0x7ffffff) << 3 # (ns)
284  for finesse in range(0, 4):
285  dc = (finesse << 2) + (copper & 0x3)
286  sectorFB = self.dcToSectorFBdcToSectorFB[dc]
287  nWords = rawklm.GetDetectorNwords(0, finesse)
288  if nWords <= 0:
289  continue
290  bufSlot = rawklm.GetDetectorBuffer(0, finesse)
291  lastWord = bufSlot[nWords - 1]
292  windowValue = (lastWord >> 16) & 0xffff
293  if windowValue != self.windowValuewindowValue:
294  if windowValue in self.dict_nEventsdict_nEvents:
295  self.windowValuewindowValue = windowValue
296  self.hist_mappedRPCTimeCalByWindowhist_mappedRPCTimeCalByWindow = self.dict_mappedRPCTimeCalByWindowdict_mappedRPCTimeCalByWindow[windowValue]
297  self.hist_occupancyXYByWindowhist_occupancyXYByWindow = self.dict_occupancyXYByWindowdict_occupancyXYByWindow[windowValue]
298  else:
299  return # skip bogus event, incuding event with seed value of 0xcafe
300  if lastWord & 0xffff != 0:
301  print("##1 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'lastWord=', hex(lastWord))
302  if (nWords % 2) == 0:
303  print("##2 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'should be odd -- skipping')
304  continue
305  n = nWords >> 1 # number of Data-Concentrator data packets
306  self.dict_nRawKLMsdict_nRawKLMs[self.windowValuewindowValue] += n
307  # first (and only) pass over this DC's hits: histogram everything
308  for j in range(0, n):
309  word0 = bufSlot[j * 2]
310  word1 = bufSlot[j * 2 + 1]
311  ctime = word0 & 0xffff
312  channel = (word0 >> 16) & 0x7f
313  axis = (word0 >> 23) & 0x01
314  lane = (word0 >> 24) & 0x1f # 1..2 for scints, 8..20 for RPCs (=readout-board slot - 7)
315  flag = (word0 >> 30) & 0x03 # 1 for RPCs, 2 for scints
316  tdc = (word1 >> 16) & 0x07ff
317  isRPC = (flag == 1)
318  electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
319  if electId in self.electIdToModuleIdelectIdToModuleId: # BKLM mapped-channel histograms
320  if isRPC:
321  tCal = int(tdc - trigCtime) & 0x03ff # in ns, range is 0..1023
322  self.hist_mappedRPCTimeCalhist_mappedRPCTimeCal.Fill(tCal)
323  self.hist_mappedRPCTimeCalByWindowhist_mappedRPCTimeCalByWindow.Fill(tCal)
324 
325  # for normalization of the hit-counter dictionaries, now that we know that this is a valid event
326 
327  self.dict_nEventsdict_nEvents[self.windowValuewindowValue] += 1
328 
329  # Process the BKLMHit2ds
330 
331  self.dict_nHit2dsdict_nHit2ds[self.windowValuewindowValue] += len(hit2ds)
332  for hit2d in hit2ds:
333  key = hit2d.getModuleID()
334  fb = (key & self.BKLM_SECTION_MASKBKLM_SECTION_MASK) >> self.BKLM_SECTION_BITBKLM_SECTION_BIT
335  x = hit2d.getGlobalPositionX()
336  y = hit2d.getGlobalPositionY()
337  if fb == 0: # backward
338  self.hist_occupancyBackwardXYhist_occupancyBackwardXY.Fill(x, y)
339  else: # forward
340  self.hist_occupancyForwardXYhist_occupancyForwardXY.Fill(x, y)
341  self.hist_occupancyXYByWindowhist_occupancyXYByWindow.Fill(x, y)
def fillDB()
Definition: bklmDB.py:17
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:56
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
int EKLM_ID
COPPER base identifier for EKLM readout.
hist_mappedRPCTimeCalByWindow
reference to the RPC-time histogram for the currevent value of the lookback window parameter
hist_occupancyXYByWindow
reference to the xy scatterplot for the currevent value of the lookback window parameter
dict_nHit2ds
dictionary of the number of BKLMHit2ds for each lookback-window value
hist_mappedRPCTimeCal
histogram of RPC TDC - trigger value
hist_occupancyBackwardXY
scatterplot of end view of backward BKLM for all BKLMHit2ds
windowMaxValue
highest observed lookback-window value
hist_occupancyForwardXY
scatterplot of end view of forward BKLM for all BKLMHit2ds
tuple BKLM_SECTION_MASK
bit mask for section [0..1]; forward is 0
dict_mappedRPCTimeCalByWindow
dictionary of histograms of RPC TDC - trigger value, keyed by lookback-window value
def __init__(self, exp, run, histName, pdfName, mode, window)
dcToSectorFB
map for data concentrator -> sectorFB
int BKLM_ID
COPPER base identifier for BKLM readout.
dict_occupancyXYByWindow
dictionary of scatterplots of end view of forward BKLM, keyed by lookback-window value
dict_nRawKLMs
dictionary of the number of RawKLM hits for each lookback-window value
windowMinValue
highest observed lookback-window value
int BKLM_SECTION_BIT
bit position for section [0..1]; forward is 0
pdfName
internal copy of the pathname of the output histogram PDF file
histogramFile
Output ROOT TFile that will contain the histograms/scatterplots.
histName
internal copy of the pathname of the output histogram ROOT file
dict_nEvents
dictionary of the number of events for each lookback-window value, for normalization
electIdToModuleId
readout <-> detector map (from the information retrieved from the conditions database)
windowMode
window mode as a string for histogram labels/titles
windowValue
cached value of the lookback-window value, to avoid unnecessary reassignments-to-same-value in event(...
sectorFBToDC
map for sectorFB -> data concentrator