Belle II Software  release-06-01-15
EventDisplayer.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # Purpose:
13 # basf module to create BKLM event displays from BKLMHit2ds, ExtHits, and MuidHits.
14 #
15 
16 import basf2
17 import bklmDB
18 import math
19 import ROOT
20 from ROOT import Belle2
21 
22 # =========================================================================
23 #
24 # EventDisplayer
25 #
26 # =========================================================================
27 
28 
29 class EventDisplayer(basf2.Module):
30  """Draw BKLM event displays from BKLMHit2ds, ExtHits, and MuidHits."""
31 
32 
33  BKLM_ID = 0x07000000
34 
35  EKLM_ID = 0x08000000
36 
37  BKLM_STRIP_BIT = 0
38 
39  BKLM_PLANE_BIT = 6
40 
41  BKLM_LAYER_BIT = 7
42 
43  BKLM_SECTOR_BIT = 11
44 
45  BKLM_SECTION_BIT = 14
46 
47  BKLM_MAXSTRIP_BIT = 15
48 
49  BKLM_STRIP_MASK = 0x3f
50 
51  BKLM_PLANE_MASK = (1 << BKLM_PLANE_BIT)
52 
53  BKLM_LAYER_MASK = (15 << BKLM_LAYER_BIT)
54 
55  BKLM_SECTOR_MASK = (7 << BKLM_SECTOR_BIT)
56 
57  BKLM_SECTION_MASK = (1 << BKLM_SECTION_BIT)
58 
59  BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
60 
61  BKLM_MODULEID_MASK = (BKLM_SECTION_MASK | BKLM_SECTOR_MASK | BKLM_LAYER_MASK)
62 
63  def __init__(self, exp, run, eventPdfName, maxDisplays, minRPCHits, minMuidHits):
64  """Constructor
65 
66  Arguments:
67  exp (str): formatted experiment number
68  run (str): formatter run number
69  eventPdfName (str): path name of the output event-display PDF file
70  maxDisplays (int): max # of events displays to write
71  minRPCHits (int): min # of RPC BKLMHit2ds in any sector for event display
72  minMuidHits (int): min # of MuidHits in the event for event display
73  """
74  super().__init__()
75 
76  self.expexp = exp
77 
78  self.runrun = run
79 
80  self.eventPdfNameeventPdfName = eventPdfName
81 
82  self.maxDisplaysmaxDisplays = maxDisplays
83 
84  self.minRPCHitsminRPCHits = minRPCHits
85 
86  self.minMuidHitsminMuidHits = minMuidHits
87 
88  self.eventCountereventCounter = 0
89 
90  self.eventDisplayseventDisplays = 0
91 
92  self.lastTitlelastTitle = ''
93 
94  def initialize(self):
95  """Handle job initialization: fill the mapping database, create histograms, open the event-display file"""
96 
97  expRun = 'e{0:02d}r{1}: '.format(int(self.expexp), int(self.runrun))
98 
99  # Open the output PDF file for event displays
100 
101 
102  self.eventCanvaseventCanvas = ROOT.TCanvas("eventCanvas", self.eventPdfNameeventPdfName, 3200, 1600)
103  title = '{0}['.format(self.eventPdfNameeventPdfName)
104  self.eventCanvaseventCanvas.SaveAs(title)
105  self.eventCanvaseventCanvas.Clear()
106  self.eventCanvaseventCanvas.Divide(2, 1)
107 
108  # Create the boilerplate for the end- and side-views of the event display
109 
110 
111  self.cosinecosine = [0, 0, 0, 0, 0, 0, 0, 0]
112 
113  self.sinesine = [0, 0, 0, 0, 0, 0, 0, 0]
114  for sector in range(0, 8):
115  phi = math.pi * sector / 4
116  self.cosinecosine[sector] = math.cos(phi)
117  self.sinesine[sector] = math.sin(phi)
118 
119  self.hist_XYhist_XY = ROOT.TH2F('XY', ' ;x;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
120  self.hist_XYhist_XY.SetStats(False)
121 
122  self.hist_ZYhist_ZY = ROOT.TH2F('ZY', ' ;z;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
123  self.hist_ZYhist_ZY.SetStats(False)
124  # 300x300 cm^2 grid for each octant
125  u1 = 65
126  u2 = 365
127  u3 = 100
128  u4 = 400
129  u5 = 150
130 
131  self.hist_XYShist_XYS = [0, 0, 0, 0, 0, 0, 0, 0]
132  self.hist_XYShist_XYS[0] = ROOT.TH2F('XYS0', ' ;x;y', 10, +u3, +u4, 10, -u5, +u5)
133  self.hist_XYShist_XYS[0].SetStats(False)
134  self.hist_XYShist_XYS[1] = ROOT.TH2F('XYS1', ' ;x;y', 10, +u1, +u2, 10, +u1, +u2)
135  self.hist_XYShist_XYS[1].SetStats(False)
136  self.hist_XYShist_XYS[2] = ROOT.TH2F('XYS2', ' ;x;y', 10, -u5, +u5, 10, +u3, +u4)
137  self.hist_XYShist_XYS[2].SetStats(False)
138  self.hist_XYShist_XYS[3] = ROOT.TH2F('XYS3', ' ;x;y', 10, -u2, -u1, 10, +u1, +u2)
139  self.hist_XYShist_XYS[3].SetStats(False)
140  self.hist_XYShist_XYS[4] = ROOT.TH2F('XYS4', ' ;x;y', 10, -u4, -u3, 10, -u5, +u5)
141  self.hist_XYShist_XYS[4].SetStats(False)
142  self.hist_XYShist_XYS[5] = ROOT.TH2F('XYS5', ' ;x;y', 10, -u2, -u1, 10, -u2, -u1)
143  self.hist_XYShist_XYS[5].SetStats(False)
144  self.hist_XYShist_XYS[6] = ROOT.TH2F('XYS6', ' ;x;y', 10, -u5, +u5, 10, -u4, -u3)
145  self.hist_XYShist_XYS[6].SetStats(False)
146  self.hist_XYShist_XYS[7] = ROOT.TH2F('XYS7', ' ;x;y', 10, +u1, +u2, 10, -u2, -u1)
147  self.hist_XYShist_XYS[7].SetStats(False)
148 
149  self.hist_ZYShist_ZYS = ROOT.TH2F('ZYS', ' ;z;y', 10, -150.0, 150.0, 10, 125.0, 425.0)
150  self.hist_ZYShist_ZYS.SetStats(False)
151  ROOT.gStyle.SetOptStat(10)
152 
153  self.bklmXYbklmXY = []
154  r0 = 201.9 + 0.5 * 4.4 # cm
155  dr = 9.1 # cm
156  tan0 = math.tan(math.pi / 8.0)
157  g = ROOT.TGraph()
158  g.SetPoint(0, -200.0, 0.0)
159  g.SetPoint(1, +200.0, 0.0)
160  g.SetLineColor(19)
161  g.SetLineWidth(1)
162  g.SetLineStyle(3)
163  self.bklmXYbklmXY.append(g)
164  g = ROOT.TGraph()
165  g.SetPoint(0, 0.0, -200.0)
166  g.SetPoint(1, 0.0, +200.0)
167  g.SetLineColor(19)
168  g.SetLineWidth(1)
169  g.SetLineStyle(3)
170  self.bklmXYbklmXY.append(g)
171  g = ROOT.TGraph()
172  g.SetPoint(0, -5.0, 0.0)
173  g.SetPoint(1, +5.0, 0.0)
174  g.SetPoint(2, 0.0, 0.0)
175  g.SetPoint(3, 0.0, +5.0)
176  g.SetPoint(4, 0.0, -5.0)
177  g.SetLineColor(1)
178  g.SetLineWidth(1)
179  self.bklmXYbklmXY.append(g)
180  for layer in range(0, 15):
181  r = r0 + layer * dr
182  x = r * tan0
183  g = ROOT.TGraph()
184  g.SetPoint(0, +r, -x)
185  g.SetPoint(1, +r, +x)
186  g.SetPoint(2, +x, +r)
187  g.SetPoint(3, -x, +r)
188  g.SetPoint(4, -r, +x)
189  g.SetPoint(5, -r, -x)
190  g.SetPoint(6, -x, -r)
191  g.SetPoint(7, +x, -r)
192  g.SetPoint(8, +r, -x)
193  if layer < 2:
194  g.SetLineColor(18)
195  else:
196  g.SetLineColor(17)
197  if (layer % 5) == 0:
198  g.SetLineStyle(3)
199  g.SetLineWidth(1)
200  self.bklmXYbklmXY.append(g)
201 
202  self.bklmZYbklmZY = []
203  rF = r0 + 14 * dr
204  x0 = r0 * tan0
205  z0 = 47.0 # cm
206  zL = 220.0 # cm
207  g = ROOT.TGraph()
208  g.SetPoint(0, -zL + z0 - 140.0, 0.0)
209  g.SetPoint(1, +zL + z0 + 70.0, 0.0)
210  g.SetLineColor(19)
211  g.SetLineWidth(1)
212  g.SetLineStyle(3)
213  self.bklmZYbklmZY.append(g)
214  g = ROOT.TGraph()
215  g.SetPoint(0, 0.0, -315.0)
216  g.SetPoint(1, 0.0, +340.0)
217  g.SetLineColor(19)
218  g.SetLineWidth(1)
219  g.SetLineStyle(3)
220  self.bklmZYbklmZY.append(g)
221  g = ROOT.TGraph()
222  g.SetPoint(0, -5.0, 0.0)
223  g.SetPoint(1, +5.0, 0.0)
224  g.SetPoint(2, 0.0, 0.0)
225  g.SetPoint(3, 0.0, +5.0)
226  g.SetPoint(4, 0.0, -5.0)
227  g.SetLineColor(1)
228  g.SetLineWidth(1)
229  self.bklmZYbklmZY.append(g)
230  g = ROOT.TGraph()
231  g.SetPoint(0, -zL + z0, +x0)
232  g.SetPoint(1, -zL + z0, +r0)
233  g.SetLineColor(18)
234  g.SetLineWidth(1)
235  g.SetLineStyle(3)
236  self.bklmZYbklmZY.append(g)
237  g = ROOT.TGraph()
238  g.SetPoint(0, -zL + z0, -x0)
239  g.SetPoint(1, -zL + z0, -r0)
240  g.SetLineColor(18)
241  g.SetLineWidth(1)
242  g.SetLineStyle(3)
243  self.bklmZYbklmZY.append(g)
244  g = ROOT.TGraph()
245  g.SetPoint(0, +zL + z0, +x0)
246  g.SetPoint(1, +zL + z0, +r0)
247  g.SetLineColor(18)
248  g.SetLineWidth(1)
249  g.SetLineStyle(3)
250  self.bklmZYbklmZY.append(g)
251  g = ROOT.TGraph()
252  g.SetPoint(0, +zL + z0, -x0)
253  g.SetPoint(1, +zL + z0, -r0)
254  g.SetLineColor(18)
255  g.SetLineWidth(1)
256  g.SetLineStyle(3)
257  self.bklmZYbklmZY.append(g)
258  g = ROOT.TGraph()
259  g.SetPoint(0, -zL + z0, r0)
260  g.SetPoint(1, +zL + z0, r0)
261  g.SetPoint(2, +zL + z0, rF)
262  g.SetPoint(3, -zL + z0, rF)
263  g.SetPoint(4, -zL + z0, r0)
264  g.SetLineColor(18)
265  g.SetLineWidth(1)
266  self.bklmZYbklmZY.append(g)
267  g = ROOT.TGraph()
268  g.SetPoint(0, -zL + z0, -r0)
269  g.SetPoint(1, +zL + z0, -r0)
270  g.SetPoint(2, +zL + z0, -rF)
271  g.SetPoint(3, -zL + z0, -rF)
272  g.SetPoint(4, -zL + z0, -r0)
273  g.SetLineColor(18)
274  g.SetLineWidth(1)
275  self.bklmZYbklmZY.append(g)
276  g = ROOT.TGraph()
277  g.SetPoint(0, -zL + z0, -x0)
278  g.SetPoint(1, +zL + z0, -x0)
279  g.SetPoint(2, +zL + z0, +x0)
280  g.SetPoint(3, -zL + z0, +x0)
281  g.SetPoint(4, -zL + z0, -x0)
282  g.SetLineColor(18)
283  g.SetLineWidth(1)
284  self.bklmZYbklmZY.append(g)
285 
286  self.bklmZYLbklmZYL = []
287  for layer in range(0, 15):
288  r = r0 + layer * dr
289  z0 = 47.0 # cm
290  zL = 220.0 # cm
291  g = ROOT.TGraph()
292  g.SetPoint(0, -zL + z0, r)
293  g.SetPoint(1, +zL + z0, r)
294  g.SetLineColor(19)
295  g.SetLineWidth(1)
296  self.bklmZYLbklmZYL.append(g)
297 
298  self.electIdToModuleIdelectIdToModuleId = bklmDB.fillDB()
299 
300  self.sectorFBToDCsectorFBToDC = [11, 15, 2, 6, 10, 14, 3, 7, 9, 13, 0, 4, 8, 12, 1, 5]
301 
302  self.dcToSectorFBdcToSectorFB = [10, 14, 2, 6, 11, 15, 3, 7, 12, 8, 4, 0, 13, 9, 5, 1]
303 
304  self.t0Calt0Cal = 312
305 
306  self.t0Cal2dt0Cal2d = 308
307 
308  self.ct0Calct0Cal = 455
309 
310  self.ct0Cal2dct0Cal2d = 520
311 
312  self.t0RPCt0RPC = [8, -14, -6, -14, -2, 10, 9, 13, 0, -10, -14, -20, 2, 6, 14, 11]
313 
314  self.ct0Scintct0Scint = [-1, -33, -46, -33, -2, 32, 51, 32, 0, -32, -45, -33, -4, 34, 45, 27]
315 
316  def terminate(self):
317  """Handle job termination: close event-display file"""
318 
319  pdfNameLast = '{0}]'.format(self.eventPdfNameeventPdfName)
320  self.eventCanvaseventCanvas.Print(pdfNameLast, self.lastTitlelastTitle)
321  print('Goodbye')
322 
323  def beginRun(self):
324  """Handle begin of run: print diagnostic message"""
325  EventMetaData = Belle2.PyStoreObj('EventMetaData')
326  print('beginRun', EventMetaData.getRun())
327 
328  def endRun(self):
329  """Handle end of run: print diagnostic message"""
330  EventMetaData = Belle2.PyStoreObj('EventMetaData')
331  print('endRun', EventMetaData.getRun())
332 
333  def event(self):
334  """Process one event: (optionally) draw event display"""
335 
336  super().return_value(self.eventDisplayseventDisplays < self.maxDisplaysmaxDisplays)
337 
338  if self.eventDisplayseventDisplays >= self.maxDisplaysmaxDisplays:
339  return
340 
341  self.eventCountereventCounter += 1
342  EventMetaData = Belle2.PyStoreObj('EventMetaData')
343  event = EventMetaData.getEvent()
344  rawklms = Belle2.PyStoreArray('RawKLMs') # to determine if BKLMHit2d is prompt
345  hit2ds = Belle2.PyStoreArray('BKLMHit2ds')
346  exthits = Belle2.PyStoreArray('ExtHits')
347  muidhits = Belle2.PyStoreArray('MuidHits')
348 
349  rawFb = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
350  rawSector = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
351  rawLayer = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
352  rawPlane = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
353  rawStrip = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
354  rawCtime = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
355 
356  for copper in range(0, len(rawklms)):
357  rawklm = rawklms[copper]
358  nodeID = rawklm.GetNodeID(0) - self.BKLM_IDBKLM_ID
359  if nodeID >= self.EKLM_IDEKLM_ID - self.BKLM_IDBKLM_ID:
360  nodeID = nodeID - (self.EKLM_IDEKLM_ID - self.BKLM_IDBKLM_ID) + 4
361  if (nodeID < 0) or (nodeID > 4): # skip EKLM nodes
362  continue
363  for finesse in range(0, 4):
364  dc = (finesse << 2) + copper
365  n = rawklm.GetDetectorNwords(0, finesse)
366  bufSlot = rawklm.GetDetectorBuffer(0, finesse)
367  if n > 0:
368  n = n >> 1
369  for j in range(0, n):
370  word0 = bufSlot[j * 2]
371  word1 = bufSlot[j * 2 + 1]
372  ctime = word0 & 0xffff
373  channel = (word0 >> 16) & 0x7f
374  axis = (word0 >> 23) & 0x01
375  lane = (word0 >> 24) & 0x1f # crate's slot number
376  flag = (word0 >> 30) & 0x03
377  fb = -1
378  sector = -1
379  layer = -1
380  plane = -1
381  strip = -1
382  electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
383  if electId in self.electIdToModuleIdelectIdToModuleId:
384  moduleId = self.electIdToModuleIdelectIdToModuleId[electId]
385  fb = (moduleId & self.BKLM_SECTION_MASKBKLM_SECTION_MASK) >> self.BKLM_SECTION_BITBKLM_SECTION_BIT
386  sector = (moduleId & self.BKLM_SECTOR_MASKBKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BITBKLM_SECTOR_BIT
387  layer = (moduleId & self.BKLM_LAYER_MASKBKLM_LAYER_MASK) >> self.BKLM_LAYER_BITBKLM_LAYER_BIT
388  plane = (moduleId & self.BKLM_PLANE_MASKBKLM_PLANE_MASK) >> self.BKLM_PLANE_BITBKLM_PLANE_BIT
389  strip = (moduleId & self.BKLM_STRIP_MASKBKLM_STRIP_MASK) >> self.BKLM_STRIP_BITBKLM_STRIP_BIT
390  rawFb[dc].append(fb)
391  rawSector[dc].append(sector)
392  rawLayer[dc].append(layer)
393  rawPlane[dc].append(plane)
394  rawStrip[dc].append(strip)
395  rawCtime[dc].append(ctime)
396 
397  tCal2d = []
398  for hit2d in hit2ds:
399  tCal2d.append(hit2d.getTime())
400 
401  # Process the ExtHits for event display
402 
403  xold = 0
404  yold = 0
405  zold = 0
406  sold = 0
407  epositions = []
408  sumx = 0
409  sumy = 0
410  sumz = 0
411  sumn = 0
412  for exthit in exthits:
413  sector = (exthit.getCopyID() & self.BKLM_SECTOR_MASKBKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BITBKLM_SECTOR_BIT
414  extPosition = exthit.getPosition()
415  x = extPosition[0]
416  y = extPosition[1]
417  z = extPosition[2]
418  dx = x - xold
419  dy = y - yold
420  dz = z - zold
421  if (dx * dx + dy * dy + dz * dz > 36) and (sumn > 0):
422  eposition = [sumx / sumn, sumy / sumn, sumz / sumn, sold]
423  epositions.append(eposition)
424  sumx = 0
425  sumy = 0
426  sumz = 0
427  sumn = 0
428  else:
429  if sumn == 0:
430  xold = x
431  yold = y
432  zold = z
433  sold = sector
434  sumx = sumx + x
435  sumy = sumy + y
436  sumz = sumz + z
437  sumn = sumn + 1
438  if sumn > 0:
439  eposition = [sumx / sumn, sumy / sumn, sumz / sumn, sold]
440  epositions.append(eposition)
441  extXYGraph = ROOT.TGraph()
442  extXYGraph.SetMarkerColor(30)
443  extXYGraph.SetMarkerSize(2.25)
444  extXYGraph.SetMarkerStyle(21)
445  extZYGraph = ROOT.TGraph()
446  extZYGraph.SetMarkerColor(30)
447  extZYGraph.SetMarkerSize(2.25)
448  extZYGraph.SetMarkerStyle(21)
449  extZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
450  for sector in range(0, 8):
451  extZYSGraph[sector] = ROOT.TGraph()
452  extZYSGraph[sector].SetMarkerColor(30)
453  extZYSGraph[sector].SetMarkerSize(2.25)
454  extZYSGraph[sector].SetMarkerStyle(21)
455  j = -1
456  for eposition in epositions:
457  j = j + 1
458  x = eposition[0]
459  y = eposition[1]
460  z = eposition[2]
461  extXYGraph.SetPoint(j, x, y)
462  extZYGraph.SetPoint(j, z, y)
463  sector = int(eposition[3])
464  nPoint = extZYSGraph[sector].GetN()
465  extZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosinecosine[sector] + y * self.sinesine[sector]))
466 
467  # Process the MuidHits for event display
468 
469  zMuids = [0, 0, 0, 0, 0, 0, 0, 0]
470  nMuids = [0, 0, 0, 0, 0, 0, 0, 0]
471  muidXYGraph = ROOT.TGraph()
472  muidXYGraph.SetMarkerColor(5)
473  muidXYGraph.SetMarkerSize(2.0)
474  muidXYGraph.SetMarkerStyle(20)
475  muidZYGraph = ROOT.TGraph()
476  muidZYGraph.SetMarkerColor(5)
477  muidZYGraph.SetMarkerSize(2.0)
478  muidZYGraph.SetMarkerStyle(20)
479  muidZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
480  for sector in range(0, 8):
481  muidZYSGraph[sector] = ROOT.TGraph()
482  muidZYSGraph[sector].SetMarkerColor(5)
483  muidZYSGraph[sector].SetMarkerSize(2.0)
484  muidZYSGraph[sector].SetMarkerStyle(20)
485  j = -1
486  for muidhit in muidhits:
487  j += 1
488  muidPosition = muidhit.getExtPosition()
489  x = muidPosition[0]
490  y = muidPosition[1]
491  z = muidPosition[2]
492  muidXYGraph.SetPoint(j, x, y)
493  muidZYGraph.SetPoint(j, z, y)
494  sector = muidhit.getSector()
495  nPoint = muidZYSGraph[sector].GetN()
496  muidZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosinecosine[sector] + y * self.sinesine[sector]))
497  if nMuids[sector] == 0:
498  zMuids[sector] = z
499  nMuids[sector] += 1
500 
501  # Process the BKLMHit2ds for event display
502 
503  rpcHitCount = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
504  promptXYGraph = ROOT.TGraph()
505  promptXYGraph.SetMarkerColor(4)
506  promptXYGraph.SetMarkerSize(2.0)
507  promptXYGraph.SetMarkerStyle(29)
508  promptZYGraph = ROOT.TGraph()
509  promptZYGraph.SetMarkerColor(4)
510  promptZYGraph.SetMarkerSize(2.0)
511  promptZYGraph.SetMarkerStyle(29)
512  bkgdXYGraph = ROOT.TGraph()
513  bkgdXYGraph.SetMarkerColor(2)
514  bkgdXYGraph.SetMarkerSize(2.0)
515  bkgdXYGraph.SetMarkerStyle(29)
516  bkgdZYGraph = ROOT.TGraph()
517  bkgdZYGraph.SetMarkerColor(2)
518  bkgdZYGraph.SetMarkerSize(2.0)
519  bkgdZYGraph.SetMarkerStyle(29)
520  promptZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
521  bkgdZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
522  for sector in range(0, 8):
523  promptZYSGraph[sector] = ROOT.TGraph()
524  promptZYSGraph[sector].SetMarkerColor(4)
525  promptZYSGraph[sector].SetMarkerSize(2.0)
526  promptZYSGraph[sector].SetMarkerStyle(29)
527  bkgdZYSGraph[sector] = ROOT.TGraph()
528  bkgdZYSGraph[sector].SetMarkerColor(2)
529  bkgdZYSGraph[sector].SetMarkerSize(2.0)
530  bkgdZYSGraph[sector].SetMarkerStyle(29)
531  jPrompt = -1
532  jBkgd = -1
533  for hit2d in hit2ds:
534  key = hit2d.getModuleID()
535  layer = (key & self.BKLM_LAYER_MASKBKLM_LAYER_MASK) >> self.BKLM_LAYER_BITBKLM_LAYER_BIT
536  sector = (key & self.BKLM_SECTOR_MASKBKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BITBKLM_SECTOR_BIT
537  fb = (key & self.BKLM_SECTION_MASKBKLM_SECTION_MASK) >> self.BKLM_SECTION_BITBKLM_SECTION_BIT
538  phiStripMin = hit2d.getPhiStripMin() - 1
539  phiStripMax = hit2d.getPhiStripMax() - 1
540  zStripMin = hit2d.getZStripMin() - 1
541  zStripMax = hit2d.getZStripMax() - 1
542  sectorFB = sector if fb == 0 else sector + 8
543  if layer >= 2:
544  rpcHitCount[sectorFB] += 1
545  dc = self.sectorFBToDCsectorFBToDC[sectorFB]
546  copper = dc & 0x03
547  finesse = dc >> 2
548  n = rawklms[copper].GetDetectorNwords(0, finesse) >> 1
549  trigCtime = (rawklms[copper].GetTTCtime(0) & 0x07ffffff) << 3
550  ctDiffMax = 99999
551  tCal = -1
552  jZ = -1
553  jPhi = -1
554  ctZ = 0
555  ctPhi = 0
556  for j in range(0, n):
557  if layer != rawLayer[dc][j]:
558  continue
559  if sector != rawSector[dc][j]:
560  continue
561  if fb != rawFb[dc][j]:
562  continue
563  strip = rawStrip[dc][j]
564  plane = rawPlane[dc][j]
565  if plane == 0: # it's a z strip
566  if strip < zStripMin:
567  continue
568  if strip > zStripMax:
569  continue
570  ctZ = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
571  jZ = j
572  else: # it's a phi strip
573  if strip < phiStripMin:
574  continue
575  if strip > phiStripMax:
576  continue
577  ctPhi = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
578  jPhi = j
579  if (jZ >= 0) and (jPhi >= 0):
580  if layer < 2: # it's a scint layer
581  if abs(ctZ - ctPhi) > 40:
582  continue
583  ct = int((ctZ + ctPhi) * 0.5 - trigCtime - self.ct0Scintct0Scint[sectorFB]) & 0x3ff
584  if abs(ct - self.ct0Calct0Cal) < ctDiffMax:
585  ctDiffMax = int(abs(ct - self.ct0Calct0Cal))
586  tCal = ct
587  if ctDiffMax == 0:
588  break
589  else: # it's an RPC layer
590  tCal = ((int(hit2d.getTime()) - trigCtime) & 0x03ff) - self.t0RPCt0RPC[sectorFB] - 0.75 * jPhi - 0.75 * jZ
591  break
592  x = hit2d.getGlobalPositionX()
593  y = hit2d.getGlobalPositionY()
594  z = hit2d.getGlobalPositionZ()
595  isPromptHit = False
596  if layer < 2:
597  if abs(tCal - self.ct0Cal2dct0Cal2d) < 20:
598  isPromptHit = True
599  else:
600  if abs(tCal - self.t0Cal2dt0Cal2d) < 20:
601  isPromptHit = True
602  if isPromptHit:
603  jPrompt += 1
604  promptXYGraph.SetPoint(jPrompt, x, y)
605  promptZYGraph.SetPoint(jPrompt, z, y)
606  nPoint = promptZYSGraph[sector].GetN()
607  promptZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosinecosine[sector] + y * self.sinesine[sector]))
608  else:
609  jBkgd += 1
610  bkgdXYGraph.SetPoint(jBkgd, x, y)
611  bkgdZYGraph.SetPoint(jBkgd, z, y)
612  nPoint = bkgdZYSGraph[sector].GetN()
613  bkgdZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosinecosine[sector] + y * self.sinesine[sector]))
614 
615  hasEnoughRPCHits = False
616  for count in rpcHitCount:
617  if count > self.minRPCHitsminRPCHits:
618  hasEnoughRPCHits = True
619  break
620  if hasEnoughRPCHits and (len(muidhits) > self.minMuidHitsminMuidHits):
621  self.eventDisplayseventDisplays += 1
622  title = 'e{0:02d}r{1}: event {2}'.format(int(self.expexp), int(self.runrun), event)
623  self.hist_XYhist_XY.SetTitle(title)
624  self.hist_ZYhist_ZY.SetTitle(title)
625  self.eventCanvaseventCanvas.cd(1)
626  self.hist_XYhist_XY.Draw()
627  for g in self.bklmXYbklmXY:
628  g.Draw("L")
629  if extXYGraph.GetN() > 0:
630  extXYGraph.Draw("P")
631  if muidXYGraph.GetN() > 0:
632  muidXYGraph.Draw("P")
633  if bkgdXYGraph.GetN() > 0:
634  bkgdXYGraph.Draw("P")
635  if promptXYGraph.GetN() > 0:
636  promptXYGraph.Draw("P")
637  self.eventCanvaseventCanvas.cd(2)
638  self.hist_ZYhist_ZY.Draw()
639  for g in self.bklmZYbklmZY:
640  g.Draw("L")
641  if extZYGraph.GetN() > 0:
642  extZYGraph.Draw("P")
643  if muidZYGraph.GetN() > 0:
644  muidZYGraph.Draw("P")
645  if bkgdZYGraph.GetN() > 0:
646  bkgdZYGraph.Draw("P")
647  if promptZYGraph.GetN() > 0:
648  promptZYGraph.Draw("P")
649  self.lastTitlelastTitle = "Title:E{0} (#{1})".format(event, self.eventCountereventCounter)
650  self.eventCanvaseventCanvas.Print(self.eventPdfNameeventPdfName, self.lastTitlelastTitle)
651  for sector in range(0, 8):
652  if nMuids[sector] > 0:
653  title = 'e{0:02d}r{1}: event {2} sector {3}'.format(int(self.expexp), int(self.runrun), event, sector)
654  self.hist_XYShist_XYS[sector].SetTitle(title)
655  self.eventCanvaseventCanvas.cd(1)
656  self.hist_XYShist_XYS[sector].Draw()
657  for g in self.bklmXYbklmXY:
658  g.Draw("L")
659  if extXYGraph.GetN() > 0:
660  extXYGraph.Draw("P")
661  if muidXYGraph.GetN() > 0:
662  muidXYGraph.Draw("P")
663  if bkgdXYGraph.GetN() > 0:
664  bkgdXYGraph.Draw("P")
665  if promptXYGraph.GetN() > 0:
666  promptXYGraph.Draw("P")
667  self.eventCanvaseventCanvas.cd(2)
668  z0 = zMuids[sector]
669  self.hist_ZYShist_ZYS.SetTitle(title)
670  self.hist_ZYShist_ZYS.SetBins(10, z0 - 150.0, z0 + 150.0, 10, 125.0, 425.0)
671  self.hist_ZYShist_ZYS.Draw()
672  for g in self.bklmZYbklmZY:
673  g.Draw("L")
674  for g in self.bklmZYLbklmZYL:
675  g.Draw("L")
676  if extZYSGraph[sector].GetN() > 0:
677  extZYSGraph[sector].Draw("P")
678  if muidZYSGraph[sector].GetN() > 0:
679  muidZYSGraph[sector].Draw("P")
680  if bkgdZYSGraph[sector].GetN() > 0:
681  bkgdZYSGraph[sector].Draw("P")
682  if promptZYSGraph[sector].GetN() > 0:
683  promptZYSGraph[sector].Draw("P")
684  self.lastTitlelastTitle = "Title:E{0} sector {1}".format(event, sector)
685  self.eventCanvaseventCanvas.Print(self.eventPdfNameeventPdfName, self.lastTitlelastTitle)
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.
run
internal copy of run number
hist_ZYS
blank scatterplot to define the per-sector bounds of the rotated BKLM side view
sine
table of sines for the BKLM-sector normals
hist_XY
blank scatterplot to define the bounds of the BKLM end view
int BKLM_SECTOR_BIT
bit position for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
ct0Scint
per-sector variations in scint-ctime calibration adjustment (ns) for rawKLMs
ct0Cal
scint-ctime calibration adjustment (ns) for rawKLMs
tuple BKLM_SECTOR_MASK
bit mask for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
ct0Cal2d
scint-ctime calibration adjustment (ns) for BKLMHit2ds
tuple BKLM_SECTION_MASK
bit mask for section [0..1]; forward is 0
eventCanvas
TCanvas on which event displays will be drawn.
def __init__(self, exp, run, eventPdfName, maxDisplays, minRPCHits, minMuidHits)
eventPdfName
internal copy of the pathname of the output event-display PDF file
maxDisplays
internal copy of the maximum number of event displays to write
t0RPC
per-sector variations in RPC-time calibration adjustment (ns) for rawKLMs
int BKLM_PLANE_BIT
bit position for plane-1 [0..1]; 0 is inner-plane
eventCounter
event counter (needed for PDF table of contents' ordinal event#)
bklmZYL
list of line-segment (z,y) points for the BKLM sector's zoomed and rotated side view
t0Cal
RPC-time calibration adjustment (ns) for rawKLMs.
minMuidHits
internal copy of the minimum number of MuidHits in the event for event display
dcToSectorFB
map for data concentrator -> sectorFB
int BKLM_STRIP_BIT
bit position for strip-1 [0..47]
hist_XYS
list of blank scatterplots to define the per-sector bounds of the BKLM end view
int BKLM_ID
COPPER base identifier for BKLM readout.
hist_ZY
blank scatterplot to define the bounds of the BKLM side view
minRPCHits
internal copy of the minimum number of RPC BKLMHit2ds in any sector for event display
exp
internal copy of experiment number
tuple BKLM_LAYER_MASK
bit mask for layer-1 [0..15]; 0 is innermost and 14 is outermost
bklmZY
list of line-segment (z,y) points for the BKLM side view
cosine
table of cosines for the BKLM-sector normals
eventDisplays
event-display counter
lastTitle
title of the last-drawn event display (needed for PDF table of contents' last event)
tuple BKLM_PLANE_MASK
bit mask for plane-1 [0..1]; 0 is inner-plane
int BKLM_SECTION_BIT
bit position for section [0..1]; forward is 0
int BKLM_LAYER_BIT
bit position for layer-1 [0..14]; 0 is innermost
bklmXY
list of line-segment (x,y) points for the BKLM end view
t0Cal2d
RPC-time calibration adjustment (ns) for BKLMHit2ds.
electIdToModuleId
readout <-> detector map (from the information retrieved from the conditions database)
int BKLM_STRIP_MASK
bit mask for strip-1 [0..47]
sectorFBToDC
map for sectorFB -> data concentrator