Belle II Software  release-05-01-25
EventInspector.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 # Purpose:
5 # basf module to histogram useful values in RawKLM, KLMDigit, BKLMHit1d, and BKLMHit2d
6 # data-objects in a DST ROOT file and to create BKLM event displays from these data-objects.
7 #
8 
9 import basf2
10 import bklmDB
11 import math
12 import ctypes
13 import ROOT
14 from ROOT import Belle2, TH1F, TH2F, TCanvas, THistPainter, TPad
15 
16 
17 class EventInspector(basf2.Module):
18  """Fill BKLM histograms of values from RawKLMs, KLMDigits, BKLMHit1ds, and BKLMHit2ds;
19  (optionally) draw event displays from these data-objects."""
20 
21 
22  BKLM_ID = 0x07000000
23 
24  EKLM_ID = 0x08000000
25 
26  BKLM_STRIP_BIT = 0
27 
28  BKLM_PLANE_BIT = 6
29 
30  BKLM_LAYER_BIT = 7
31 
32  BKLM_SECTOR_BIT = 11
33 
34  BKLM_SECTION_BIT = 14
35 
36  BKLM_MAXSTRIP_BIT = 15
37 
38  BKLM_STRIP_MASK = 0x3f
39 
40  BKLM_PLANE_MASK = (1 << BKLM_PLANE_BIT)
41 
42  BKLM_LAYER_MASK = (15 << BKLM_LAYER_BIT)
43 
44  BKLM_SECTOR_MASK = (7 << BKLM_SECTOR_BIT)
45 
46  BKLM_SECTION_MASK = (1 << BKLM_SECTION_BIT)
47 
48  BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
49 
50  BKLM_MODULEID_MASK = (BKLM_SECTION_MASK | BKLM_SECTOR_MASK | BKLM_LAYER_MASK)
51 
52  def __init__(self, exp, run, histName, pdfName, eventPdfName, verbosity,
53  maxDisplays, minRPCHits, legacyTimes, singleEntry, view):
54  """Constructor
55 
56  Arguments:
57  exp (str): formatted experiment number
58  run (str): formatter run number
59  histName (str): path name of the output histogram ROOT file
60  pdfName (str): path name of the output histogram PDF file
61  eventPdfName (str): path name of the output event-display PDF file
62  verbosity (int): determines how many histograms are written to the histogram PDF file
63  maxDisplays (int): max # of events displays to write
64  minRPCHits (int): min # of RPC BKLMHit2ds in any sector for event display
65  legacyTimes (bool): true to correct BKLMHit{1,2}d times in legacy reconstruction, False otherwise
66  singleEntry (int): select events with any (0) or exactly one (1) or more than one (2) entries/channel
67  view (int): view event displays using one-dimensional (1) or two-dimensional (2) BKLMHits
68  """
69  super().__init__()
70 
71  self.exp = exp
72 
73  self.run = run
74 
75  self.histName = histName
76 
77  self.pdfName = pdfName
78 
79  self.eventPdfName = eventPdfName
80 
81  self.verbosity = verbosity
82 
83  self.maxDisplays = maxDisplays
84 
85  self.minRPCHits = minRPCHits
86 
87  self.legacyTimes = legacyTimes
88 
89  self.singleEntry = singleEntry
90 
91  self.view = view
92 
93  self.eventCounter = 0
94 
95  self.eventDisplays = 0
96 
97  self.lastTitle = ''
98 
99  def makeGraph(self, x, y):
100  """Create and return a ROOT TGraph
101 
102  Arguments:
103  x[] (real): x coordinates
104  y[] (real): y coordinates
105  """
106  graph = ROOT.TGraph()
107  for i in range(0, len(x)):
108  graph.SetPoint(i, x[i], y[i])
109  graph.SetLineColor(2)
110  graph.SetLineWidth(1)
111  return graph
112 
113  def makeText(self, x, y, s):
114  """Create and return a ROOT TLatex with the following properties:
115  size = 0.04, color = red, alignment = middle centre, angle = 90 degrees
116 
117  Arguments:
118  x (real): x coordinate
119  y (real): y coordinate
120  s (str): character string
121  """
122  text = ROOT.TLatex(x, y, s)
123  text.SetTextSize(0.04)
124  text.SetTextColor(2)
125  text.SetTextAlign(22)
126  text.SetTextAngle(90)
127  return text
128 
129  def initialize(self):
130  """Handle job initialization: fill the mapping database, create histograms, open the event-display file"""
131 
132  expRun = 'e{0:02d}r{1}: '.format(int(self.exp), int(self.run))
133 
134  self.hist_XY = ROOT.TH2F('XY', ' ;x;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
135  self.hist_XY.SetStats(False)
136 
137  self.hist_ZY1D = [0, 0]
138  self.hist_ZY1D[0] = ROOT.TH2F('ZY0', ' ;z;y', 10, -200.0, 300.0, 10, -150.0, 350.0)
139  self.hist_ZY1D[1] = ROOT.TH2F('ZY1', ' ;z;y', 10, -200.0, 300.0, 10, -150.0, 350.0)
140  self.hist_ZY1D[0].SetStats(False)
141  self.hist_ZY1D[0].SetStats(False)
142 
143  self.hist_ZY = ROOT.TH2F('ZY', ' ;z;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
144  self.hist_ZY.SetStats(False)
145 
146  # All histograms/scatterplots in the output file will show '# of events' only
147  ROOT.gStyle.SetOptStat(10)
148 
150 
151  self.sectorFBToDC = [11, 15, 2, 6, 10, 14, 3, 7, 9, 13, 0, 4, 8, 12, 1, 5]
152 
153  self.dcToSectorFB = [10, 14, 2, 6, 11, 15, 3, 7, 12, 8, 4, 0, 13, 9, 5, 1]
154 
156  self.t0Cal = 312
157 
158  self.t0Cal1d = 325
159 
160  self.t0Cal2d = 308
161 
162  self.ct0Cal = 455
163 
164  self.ct0Cal1d = 533
165 
166  self.ct0Cal2d = 520
167 
168  self.t0RPC = [[[0, 0, 17, 16, 15, 13, 14, 14, 16, 17, 18, 16, 16, 18, 16],
169  [0, 0, 0, -2, -1, -5, -6, -5, -5, -3, -4, -6, -8, -6, -8],
170  [0, 0, 6, 3, 5, 0, -1, -1, -2, 0, 0, -2, -3, -2, -4],
171  [0, 0, -4, -3, -6, -8, -8, -9, -9, -7, -7, -10, -10, -10, -12],
172  [0, 0, 6, 8, 4, 4, 3, 8, 4, 7, 7, 3, 5, 5, 4],
173  [0, 0, 18, 20, 18, 18, 16, 21, 17, 20, 20, 19, 21, 19, 20],
174  [0, 0, 19, 19, 19, 18, 17, 18, 22, 24, 25, 22, 22, 24, 22],
175  [0, 0, 19, 19, 18, 17, 16, 18, 18, 20, 21, 19, 20, 20, 20],
176  [0, 0, 6, 7, 9, 5, 4, 6, 6, 9, 9, 6, 7, 8, 7],
177  [0, 0, 3, 2, 2, -4, -1, -2, -2, 1, 0, -4, 246, -3, -4],
178  [0, 0, -1, -1, -1, -5, -6, -6, -8, -5, -4, -7, -7, -7, -8],
179  [0, 0, -5, -5, -5, -12, -9, -10, -8, -6, -10, -8, -8, -11, -12],
180  [0, 0, 12, 12, 13, 8, 6, 6, 7, 9, 10, 7, 6, 8, 5],
181  [0, 0, 14, 15, 43, 12, 10, 12, 11, 15, 16, 28, 14, 15, 14],
182  [0, 0, 22, 22, 21, 19, 19, 19, 21, 23, 24, 21, 22, 22, 22],
183  [0, 0, 18, 18, 17, 16, 16, 18, 18, 20, 21, 18, 18, 20, 19]],
184  [[0, 0, 6, 5, 4, 1, 1, 2, 3, 2, 2, -1, 0, 1, -1],
185  [0, 0, -11, -12, -11, -15, -18, -18, -18, -18, -19, -22, -23, -22, -24],
186  [0, 0, -4, -7, -6, -11, -12, -12, -14, -15, -15, -18, -19, -18, -20],
187  [0, 0, -15, -15, -16, -19, -22, -21, -22, -22, -22, -25, -26, -26, -27],
188  [0, 0, -5, -3, -6, -7, -9, -9, -9, -8, -8, -13, -12, -10, -13],
189  [0, 0, 6, 7, 5, 5, 3, 9, 4, 5, 6, 3, 5, 3, 4],
190  [0, 0, 9, 10, 10, 7, 7, 7, 9, 9, 9, 6, 6, 8, 8],
191  [0, 0, 7, 8, 7, 6, 4, 5, 4, 5, 5, 4, 3, 4, 3],
192  [0, 0, -5, -3, -1, -4, -8, -7, -7, -6, -6, -6, -9, -9, -9],
193  [0, 0, -8, -8, -11, -10, -14, -15, -16, -14, -15, -20, -20, -13, -20],
194  [0, 0, -12, -12, -14, -16, -16, -15, -21, -19, -19, -23, -23, -23, -24],
195  [0, 0, -15, -15, -15, -21, -22, -22, -22, -21, -23, -25, -25, -26, -27],
196  [0, 0, 0, 0, 2, -4, -5, -5, -4, -2, -1, -5, -5, -3, -7],
197  [0, 0, 3, 3, 32, 1, 0, -1, -3, 2, 1, 13, -1, 0, -2],
198  [0, 0, 11, 11, 10, 9, 6, 7, 6, 8, 8, 5, 6, 7, 6],
199  [0, 0, 7, 8, 7, 5, 3, 5, 7, 5, 5, 2, 7, 4, 3]]]
200 
201 
202  self.ct0Scint = [[[5, 7], [-27, -24], [-29, -45], [-27, -32], [3, 6], [34, 35], [48, 44], [33, 38],
203  [4, 7], [-28, -27], [-39, -34], [-36, -33], [2, 5], [25, 30], [46, 49], [41, 31]],
204  [[0, 0], [-32, -32], [-29, -54], [-31, -40], [-1, -1], [29, 27], [41, 41], [28, 28],
205  [-2, -1], [-32, -34], [-38, -45], [-40, -41], [-3, -3], [21, 20], [41, 42], [41, 19]]]
206 
207 
208  self.histogramFile = ROOT.TFile.Open(self.histName, "RECREATE")
209  # All histograms/scatterplots in the output file will show '# of events' only
210  ROOT.gStyle.SetOptStat(10)
211  ROOT.gStyle.SetOptFit(111)
212 
213  # create the rawKLM histograms
214 
215 
216  self.hist_nDigit = ROOT.TH1F('NDigit', expRun + '# of BKLMDigits', 500, -0.5, 499.5)
217 
218  self.hist_nRawKLM = ROOT.TH1F('NRawKLM', expRun + '# of RawKLMs', 10, -0.5, 9.5)
219 
220  self.hist_rawKLMnumEvents = ROOT.TH1F('RawKLMnumEvents', expRun + 'RawKLM NumEvents;(should be 1)', 10, -0.5, 9.5)
221 
222  self.hist_rawKLMnumNodes = ROOT.TH1F('RawKLMnumNodes', expRun + 'RawKLM NumNodes;(should be 1)', 10, -0.5, 9.5)
223 
224  self.hist_rawKLMnodeID = ROOT.TH2F('RawKLMnodeID',
225  expRun + 'RawKLM NodeID;' +
226  'NodeID (bklm: 1..4, eklm:5..8);' +
227  'Copper index',
228  10, -0.5, 9.5, 10, -0.5, 9.5)
229 
230  self.hist_rawKLMlaneFlag = ROOT.TH2F('rawKLMlaneFlag',
231  expRun + 'RawKLM lane vs flag;' +
232  'Flag (1=RPC, 2=Scint);' +
233  'Lane (scint: 1..7, RPC: 8..20)',
234  4, -0.5, 3.5, 21, -0.5, 20.5)
235 
236  self.hist_rawKLMtdcExtraRPC = ROOT.TH2F('rawKLMtdcExtraRPC',
237  expRun + 'RawKLM RPC tdcExtra bits;' +
238  'Sector # (0-7 = backward, 8-15 = forward);' +
239  'tdcExtra [should be 0]',
240  16, -0.5, 15.5, 32, -0.5, 31.5)
241 
242  self.hist_rawKLMadcExtraRPC = ROOT.TH2F('rawKLMadcExtraRPC',
243  expRun + 'RawKLM RPC adcExtra bits;' +
244  'Sector # (0-7 = backward, 8-15 = forward);' +
245  'adcExtra [should be 0]',
246  16, -0.5, 15.5, 16, -0.5, 15.5)
247 
248  self.hist_rawKLMtdcExtraScint = ROOT.TH2F('rawKLMtdcExtraScint',
249  expRun + 'RawKLM Scint tdcExtra bits;' +
250  'Sector # (0-7 = backward, 8-15 = forward);' +
251  'tdcExtra',
252  16, -0.5, 15.5, 32, -0.5, 31.5)
253 
254  self.hist_rawKLMadcExtraScint = ROOT.TH2F('rawKLMadcExtraScint',
255  expRun + 'RawKLM Scint adcExtra bits;' +
256  'Sector # (0-7 = backward, 8-15 = forward);' +
257  'adcExtra',
258  16, -0.5, 15.5, 16, -0.5, 15.5)
259 
260  self.hist_rawKLMsizeMultihit = ROOT.TH1F('rawKLMsizeMultihit', expRun + 'RawKLM word count (N/channel)', 400, -0.5, 799.5)
261 
262  self.hist_rawKLMsize = ROOT.TH1F('rawKLMsize', expRun + 'RawKLM word count (1/channel)', 250, -0.5, 499.5)
263 
265 
267  for sectorFB in range(0, 16):
268  dc = self.sectorFBToDC[sectorFB]
269  copper = dc & 0x03
270  finesse = dc >> 2
271  label = 'rawKLM_S{0:02d}_sizeMultihit'.format(sectorFB)
272  title = '{0}sector {1} [COPPER {2} finesse {3}] word count (N/channel)'.format(expRun, sectorFB, copper, finesse)
273  self.hist_rawKLMsizeByDCMultihit.append(ROOT.TH1F(label, title, 100, -0.5, 199.5))
274  label = 'rawKLM_S{0:02d}_size'.format(sectorFB)
275  title = '{0}sector {1} [COPPER {2} finesse {3}] word count (1/channel)'.format(expRun, sectorFB, copper, finesse)
276  self.hist_rawKLMsizeByDC.append(ROOT.TH1F(label, title, 100, -0.5, 199.5))
277 
279 
281  for sectorFB in range(0, 16):
282  dc = self.sectorFBToDC[sectorFB]
283  copper = dc & 0x03
284  finesse = dc >> 2
285  label = 'rawKLM_S{0:02d}_channelMultiplicity'.format(sectorFB)
286  title = '{0}sector {1} [COPPER {2} finesse {3}] per-channel multiplicity (N/channel > 1);'.format(
287  expRun, sectorFB, copper, finesse) + 'Per-channel multiplicity;(Lane #) * 2 + (Axis #)'
288  self.hist_rawKLMchannelMultiplicity.append(ROOT.TH2F(label, title, 30, -0.5, 29.5, 42, -0.5, 41.5))
289  label = 'rawKLM_S{0:02d}_channelMultiplicityFine'.format(sectorFB)
290  title = '{0}sector {1} [COPPER {2} finesse {3}] per-channel multiplicity (N/channel > 1);'.format(
291  expRun, sectorFB, copper, finesse) + 'Per-channel multiplicity;(Lane #) * 256 + (Axis #) * 128 + (Channel #)'
292  self.hist_rawKLMchannelMultiplicityFine.append(ROOT.TH2F(label, title, 30, -0.5, 29.5, 8192, -0.5, 8191.5))
293 
295  'mappedSectorOccupancyMultihit',
296  expRun + 'Sector occupancy of mapped channels (N/channel);' +
297  'Sector # (0-7 = backward, 8-15 = forward)',
298  16, -0.5, 15.5)
299 
301  'unmappedSectorOccupancyMultihit',
302  expRun + 'Sector occupancy of unmapped channels (N/channel);' +
303  'Sector # (0-7 = backward, 8-15 = forward)',
304  16, -0.5, 15.5)
305 
306  self.hist_mappedSectorOccupancy = ROOT.TH1F(
307  'mappedSectorOccupancy',
308  expRun + 'Sector occupancy of mapped channels (1/channel);' +
309  'Sector # (0-7 = backward, 8-15 = forward)',
310  16, -0.5, 15.5)
311 
313  'unmappedSectorOccupancy',
314  expRun + 'Sector occupancy of unmapped channels (1/channel);' +
315  'Sector # (0-7 = backward, 8-15 = forward)',
316  16, -0.5, 15.5)
317 
319  'mappedRPCSectorOccupancy',
320  expRun + 'Sector occupancy of mapped RPC channels (1/channel);' +
321  'Sector # (0-7 = backward, 8-15 = forward)',
322  16, -0.5, 15.5)
323 
325  'mappedRPCLaneAxisOccupancy',
326  expRun + 'Lane/axis occupancy of mapped RPC channels (1/channel);' +
327  'Sector # (0-7 = backward, 8-15 = forward);' +
328  '(Lane #) * 2 + (Axis #)',
329  16, -0.5, 15.5, 42, -0.5, 41.5)
330 
332  'unmappedRPCSectorOccupancy',
333  expRun + 'Sector occupancy of unmapped RPC channels (1/channel);' +
334  'Sector # (0-7 = backward, 8-15 = forward)',
335  16, -0.5, 15.5)
336 
338  'unmappedRPCLaneAxisOccupancy',
339  expRun + 'Lane/axis occupancy of unmapped RPC channels (1/channel);' +
340  'Sector # (0-7 = backward, 8-15 = forward);' +
341  '(Lane #) * 2 + (Axis #)',
342  16, -0.5, 15.5, 42, -0.5, 41.5)
343 
345  'mappedScintSectorOccupancy',
346  expRun + 'Sector occupancy of mapped scint channels (1/channel);' +
347  'Sector # (0-7 = backward, 8-15 = forward)',
348  16, -0.5, 15.5)
349 
351  'mappedScintLaneAxisOccupancy',
352  expRun + 'Lane/axis occupancy of mapped scint channels (1/channel);' +
353  'Sector # (0-7 = backward, 8-15 = forward);' +
354  '(Lane #) * 2 + (Axis #)',
355  16, -0.5, 15.5, 42, -0.5, 41.5)
356 
358  'unmappedScintSectorOccupancy',
359  expRun + 'Sector occupancy of unmapped scint channels (1/channel);' +
360  'Sector # (0-7 = backward, 8-15 = forward)',
361  16, -0.5, 15.5)
362 
364  'unmappedScintLaneAxisOccupancy',
365  expRun + 'Lane/axis occupancy of unmapped scint channels (1/channel);' +
366  'Sector # (0-7 = backward, 8-15 = forward);' +
367  '(Lane #) * 2 + (Axis #)',
368  16, -0.5, 15.5, 42, -0.5, 41.5)
369 
371  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0],
372  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0]]
373 
375  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0],
376  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0]]
377 
379  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0],
380  [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0], [0, 0]]
381  for sectorFB in range(0, 16):
382  label = 'mappedChannelOccupancy_S{0:02d}ZPrompt'.format(sectorFB)
383  title = '{0}In-time mapped channel occupancy for sector {1} z hits;lane;channel'.format(expRun, sectorFB)
384  self.hist_mappedChannelOccupancyPrompt[sectorFB][0] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
385  label = 'mappedChannelOccupancy_S{0:02d}ZBkgd'.format(sectorFB)
386  title = '{0}Out-of-time mapped channel occupancy for sector {1} z hits;lane;channel'.format(expRun, sectorFB)
387  self.hist_mappedChannelOccupancyBkgd[sectorFB][0] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
388  label = 'unmappedChannelOccupancy_S{0:02d}Z'.format(sectorFB)
389  title = '{0}Unmapped channel occupancy for sector {1} z hits;lane;channel'.format(expRun, sectorFB)
390  self.hist_unmappedChannelOccupancy[sectorFB][0] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
391  label = 'mappedChannelOccupancy_S{0:02d}PhiPrompt'.format(sectorFB)
392  title = '{0}In-time mapped occupancy for sector {1} phi hits;lane;channel'.format(expRun, sectorFB)
393  self.hist_mappedChannelOccupancyPrompt[sectorFB][1] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
394  label = 'mappedChannelOccupancy_S{0:02d}PhiBkgd'.format(sectorFB)
395  title = '{0}Out-of-time mapped occupancy for sector {1} phi hits;lane;channel'.format(expRun, sectorFB)
396  self.hist_mappedChannelOccupancyBkgd[sectorFB][1] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
397  label = 'unmappedChannelOccupancy_S{0:02d}Phi'.format(sectorFB)
398  title = '{0}Unmapped channel occupancy for sector {1} phi hits;lane;channel'.format(expRun, sectorFB)
399  self.hist_unmappedChannelOccupancy[sectorFB][1] = ROOT.TH2F(label, title, 42, -0.25, 20.75, 128, -0.25, 63.75)
400 
401  self.hist_RPCTimeLowBitsBySector = ROOT.TH2F('RPCTimeLowBitsBySector',
402  expRun + 'RPC TDC lowest-order bits;' +
403  'Sector # (0-7 = backward, 8-15 = forward);' +
404  'TDC % 4 (ns) [should be 0]',
405  16, -0.5, 15.5, 4, -0.5, 3.5)
406 
407  self.hist_mappedRPCTime = ROOT.TH1F(
408  'mappedRPCTime', expRun + 'RPC mapped-strip time distribution;t - t(trigger) (ns)', 256, -0.5, 1023.5)
409 
410  self.hist_mappedRPCTimeCal = ROOT.TH1F(
411  'mappedRPCTimeCal', expRun + 'RPC mapped-strip time distribution;t - t(trigger) - dt(layer) (ns)', 256, -0.5, 1023.5)
412 
414 
416  for sectorFB in range(0, 16):
417  self.hist_mappedRPCPhiTimePerLayer.append([])
418  self.hist_mappedRPCZTimePerLayer.append([])
419  for layer in range(0, 15):
420  label = 'mappedRPCPhiTime_S{0:02d}L{1:02d}'.format(sectorFB, layer)
421  title = '{0}RPC sector {1} layer {2} phi time distribution;t - t(trigger) (ns)'.format(expRun, sectorFB, layer)
422  self.hist_mappedRPCPhiTimePerLayer[sectorFB].append(ROOT.TH1F(label, title, 256, -0.5, 1023.5))
423  label = 'mappedRPCZTime_S{0:02d}L{1:02d}'.format(sectorFB, layer)
424  title = '{0}RPC sector {1} layer {2} z time distribution;t - t(trigger) (ns)'.format(expRun, sectorFB, layer)
425  self.hist_mappedRPCZTimePerLayer[sectorFB].append(ROOT.TH1F(label, title, 256, -0.5, 1023.5))
426 
427  self.hist_mappedRPCTimeBySector = ROOT.TH2F('mappedRPCTimeBySector',
428  expRun + 'RPC mapped-strip time;' +
429  'Sector # (0-7 = backward, 8-15 = forward);' +
430  't - t(trigger) (ns)',
431  16, -0.5, 15.5, 128, -0.5, 1023.5)
432 
433  self.hist_mappedRPCTimeCalBySector = ROOT.TH2F('mappedRPCTimeCalBySector',
434  expRun + 'RPC mapped-strip time;' +
435  'Sector # (0-7 = backward, 8-15 = forward);' +
436  't - t(trigger) - dt(layer) (ns)',
437  16, -0.5, 15.5, 128, -0.5, 1023.5)
438 
439  self.hist_mappedRPCCtimeRange = ROOT.TH1F('mappedRPCCtimeRange',
440  expRun + 'RPC Ctime-range in event;' +
441  'CtimeMax - CtimeMin (ns)',
442  128, -0.5, 8191.5)
443 
444  self.hist_mappedRPCCtimeRangeBySector = ROOT.TH2F('mappedRPCCtimeRangeBySector',
445  expRun + 'RPC Ctime-range in event;' +
446  'Sector # (0-7 = backward, 8-15 = forward);' +
447  'CtimeMax - CtimeMin (ns)',
448  16, -0.5, 15.5, 128, -0.5, 8191.5)
449 
450  self.hist_unmappedRPCTime = ROOT.TH1F('unmappedRPCTime',
451  expRun + 'RPC unmapped-strip time distribution;' +
452  't - t(trigger) (ns)',
453  256, -0.5, 1023.5)
454 
455  self.hist_unmappedRPCTimeBySector = ROOT.TH2F('unmappedRPCTimeBySector',
456  expRun + 'RPC unmapped-strip time;' +
457  'Sector # (0-7 = backward, 8-15 = forward);' +
458  't - t(trigger) (ns)',
459  16, -0.5, 15.5, 128, -0.5, 1023.5)
460 
461  self.hist_ScintTimeLowBitsBySector = ROOT.TH2F('ScintTimeLowBitsBySector',
462  expRun + 'Scint TDC lowest-order bits;' +
463  'Sector # (0-7 = backward, 8-15 = forward);' +
464  'TDC % 4 (ns)',
465  16, -0.5, 15.5, 4, -0.5, 3.5)
466 
467  self.hist_mappedScintCtime = ROOT.TH1F('mappedScintCtime',
468  expRun + 'Scint mapped-strip ctime distribution;' +
469  'ctime - ct(trigger) (ns)',
470  32, -0.5, 1023.5)
471 
472  self.hist_mappedScintCtimeBySector = ROOT.TH2F('mappedScintCtimeBySector',
473  expRun + 'Scint mapped-strip ctime;' +
474  'Sector # (0-7 = backward, 8-15 = forward);' +
475  'ctime - ct(trigger) (ns)',
476  16, -0.5, 15.5, 32, -0.5, 1023.5)
477 
478  self.hist_mappedScintCtimeCal = ROOT.TH1F('mappedScintCtimeCal',
479  expRun + 'Scint mapped-strip ctime distribution;' +
480  'ctime - ct(trigger) - dt(layer) (ns)',
481  32, -0.5, 1023.5)
482 
483  self.hist_mappedScintCtimeCalBySector = ROOT.TH2F('mappedScintCtimeCalBySector',
484  expRun + 'Scint mapped-strip ctime;' +
485  'Sector # (0-7 = backward, 8-15 = forward);' +
486  'ctime - ct(trigger) - dt(layer) (ns)',
487  16, -0.5, 15.5, 32, -0.5, 1023.5)
488 
490 
492  for sectorFB in range(0, 16):
493  self.hist_mappedScintPhiCtimePerLayer.append([])
494  self.hist_mappedScintZCtimePerLayer.append([])
495  for layer in range(0, 2):
496  label = 'mappedScintPhiCtime_S{0:02d}L{1:02d}'.format(sectorFB, layer)
497  title = '{0}Scint sector {1} layer {2} phi ctime distribution;ctime - ct(trigger) (ns)'.format(expRun,
498  sectorFB, layer)
499  self.hist_mappedScintPhiCtimePerLayer[sectorFB].append(ROOT.TH1F(label, title, 32, -0.5, 1023.5))
500  label = 'mappedScintZCtime_S{0:02d}L{1:02d}'.format(sectorFB, layer)
501  title = '{0}Scint sector {1} layer {2} z ctime distribution;ctime - ct(trigger) (ns)'.format(expRun,
502  sectorFB, layer)
503  self.hist_mappedScintZCtimePerLayer[sectorFB].append(ROOT.TH1F(label, title, 32, -0.5, 1023.5))
504 
505  self.hist_mappedScintCtimeRange = ROOT.TH1F('mappedScintCtimeRange',
506  expRun + 'Scint ctime-range in event;' +
507  'ctimeMax - ctimeMin (ns)',
508  128, -0.5, 1023.5)
509 
510  self.hist_mappedScintCtimeRangeBySector = ROOT.TH2F('mappedScintCtimeRangeBySector',
511  expRun + 'Scint ctime-range in event;' +
512  'Sector # (0-7 = backward, 8-15 = forward);' +
513  'ctimeMax - ctimeMin (ns)',
514  16, -0.5, 15.5, 128, -0.5, 1023.5)
515 
516  self.hist_unmappedScintCtime = ROOT.TH1F('unmappedScintCtime',
517  expRun + 'Scint unmapped-strip ctime distribution;' +
518  'ctime - ct(trigger) (ns)',
519  32, -0.5, 1023.5)
520 
521  self.hist_unmappedScintCtimeBySector = ROOT.TH2F('unmappedScintCtimeBySector',
522  expRun + 'Scint unmapped-strip ctime;' +
523  'Sector # (0-7 = backward, 8-15 = forward);' +
524  'ctime - ct(trigger) (ns)',
525  16, -0.5, 15.5, 32, -0.5, 1023.5)
526 
527  self.hist_mappedScintTDC = ROOT.TH1F('mappedScintTDC',
528  expRun + 'Scint mapped-strip TDC distribution;' +
529  't (ns)',
530  32, -0.5, 31.5)
531 
532  self.hist_mappedScintTime = ROOT.TH1F('mappedScintTime',
533  expRun + 'Scint mapped-strip time distribution;' +
534  't - t(trigger) (ns)', 32, -0.5, 31.5)
535 
536  self.hist_mappedScintTDCBySector = ROOT.TH2F('mappedScintTDCBySector',
537  expRun + 'Scint mapped-strip TDC;' +
538  'Sector # (0-7 = backward, 8-15 = forward);' +
539  't (ns)',
540  16, -0.5, 15.5, 32, -0.5, 31.5)
541 
542  self.hist_mappedScintTimeBySector = ROOT.TH2F('mappedScintTimeBySector',
543  expRun + 'Scint mapped-strip time;' +
544  'Sector # (0-7 = backward, 8-15 = forward);' +
545  't - t(trigger) (ns)',
546  16, -0.5, 15.5, 32, -0.5, 31.5)
547 
548  self.hist_unmappedScintTime = ROOT.TH1F('unmappedScintTime',
549  expRun + 'Scint unmapped-strip time distribution;' +
550  't - t(trigger) (ns)',
551  32, -0.5, 31.5)
552 
553  self.hist_unmappedScintTimeBySector = ROOT.TH2F('unmappedScintTimeBySector',
554  expRun + 'Scint unmapped-strip time;' +
555  'Sector # (0-7 = backward, 8-15 = forward);' +
556  't - t(trigger) (ns)',
557  16, -0.5, 15.5, 32, -0.5, 31.5)
558 
559  # Create the RPC time-calibration/diagnostic histograms
560 
561 
562  self.hist_trigCtimeVsTrigRevo9time = ROOT.TH1F('trigCtimeVsTrigRevo9time',
563  expRun + 'trigCtime - trigRevo9time (ns)',
564  256, -1024.5, 1023.5)
565 
566  self.hist_tdcRangeRPC = ROOT.TH1F('tdcRangeRPC',
567  expRun + 'RPC TDC range;' +
568  'maxTDC - minTDC (ns)',
569  128, -0.5, 1023.5)
570 
571  self.hist_ctimeRangeRPC = ROOT.TH1F('ctimeRangeRPC',
572  expRun + 'RPC Ctime range;' +
573  'maxCtime - minCtime (ns)',
574  128, -0.5, 1023.5)
575 
576  self.hist_tdcRangeVsCtimeRangeRPC = ROOT.TH2F('tdcRangeVsCtimeRangeRPC',
577  expRun + 'RPC Ctime range vs TDC range;' +
578  'maxTDC - minTDC (ns);' +
579  'maxCtime - minCtime (ns)',
580  128, -0.5, 1023.5, 128, -0.5, 1023.5)
581 
582  self.hist_tdcRangeVsTimeRPC = ROOT.TH2F('tdcRangeVsTimeRPC',
583  expRun + 'RPC TDC range vs time;' +
584  't - t(trigger) (ns);' +
585  'maxTDC - minTDC (ns)',
586  128, -0.5, 1023.5, 128, -0.5, 1023.5)
587 
588  self.hist_ctimeRangeVsTimeRPC = ROOT.TH2F('ctimeRangeVsTimeRPC',
589  expRun + 'RPC Ctime range vs time;' +
590  't - t(trigger) (ns);' +
591  'maxCtime - minCtime (ns)',
592  128, -0.5, 1023.5, 128, -0.5, 1023.5)
593 
594  # Create the BKLMHit1d-related histograms
595 
596 
597  self.hist_nHit1d = ROOT.TH1F('NHit1d', expRun + '# of BKLMHit1ds', 100, -0.5, 99.5)
598 
599  self.hist_nHit1dRPCPrompt = ROOT.TH1F('NHit1dRPCPrompt', expRun + '# of prompt RPC BKLMHit1ds', 100, -0.5, 99.5)
600 
601  self.hist_nHit1dRPCBkgd = ROOT.TH1F('NHit1dRPCBkgd', expRun + '# of background RPC BKLMHit1ds', 100, -0.5, 99.5)
602 
603  self.hist_nHit1dScint = ROOT.TH1F('NHit1dScint', expRun + '# of scintillator BKLMHit1ds', 100, -0.5, 99.5)
604 
605  self.hist_nHit1dPrompt = ROOT.TH1F('NHit1dPrompt', expRun + '# of prompt BKLMHit1ds', 100, -0.5, 99.5)
606 
607  self.hist_nHit1dBkgd = ROOT.TH1F('NHit1dBkgd', expRun + '# of bkgd BKLMHit1ds', 100, -0.5, 99.5)
608 
609  self.hist_n1dPhiZ = ROOT.TH2F('NHit1dPhiZ',
610  expRun + 'Distribution of BKLMHit1ds;# of phi BKLMHit1ds;# of z BKLMHit1ds',
611  60, -0.5, 59.5, 60, -0.5, 59.5)
612 
613  self.hist_multiplicityPhiBySector = ROOT.TH2F('Hit1dMultiplicityPhiBySector',
614  expRun + 'BKLMHit1d phi-strip multiplicity;' +
615  'sector # (0-7 = backward, 8-15 = forward);' +
616  '# of strips',
617  16, -0.5, 15.5, 8, -0.5, 7.5)
618 
619  self.hist_multiplicityZBySector = ROOT.TH2F('Hit1dMultiplicityZBySector',
620  expRun + 'BKLMHit1d z-strip multiplicity;' +
621  'sector # (0-7 = backward, 8-15 = forward);' +
622  '# of strips',
623  16, -0.5, 15.5, 8, -0.5, 7.5)
624 
625  self.hist_tphiRPCCal1d = ROOT.TH1F('tphiRPCCal1d',
626  expRun + 'RPC BKLMHit1d phi-strip time distribution;' +
627  't(phi1D) - dt(layer) (ns)',
628  256, -0.5, 1023.5)
629 
630  self.hist_tzRPCCal1d = ROOT.TH1F('tzRPCCal1d',
631  expRun + 'RPC BKLMHit1d z-strip time distribution;' +
632  't(z1D) - dt(layer) (ns)',
633  256, -0.5, 1023.5)
634 
635  self.hist_tRPCCal1d = ROOT.TH1F('tRPCCal1d',
636  expRun + 'RPC BKLMHit1d x 2 calibrated average-time distribution;' +
637  '0.5*[t(phi1D) + t(z1D)] - dt(layer) (ns)',
638  256, -0.5, 1023.5)
639 
640  self.hist_dtRPC1d = ROOT.TH1F('dtRPC1d',
641  expRun + 'RPC BKLMHit1d x 2 time-difference distribution;' +
642  't(phi1D) - t(z1D) (ns)',
643  50, -100.0, 100.0)
644 
645  self.hist_ctphiScintCal1d = ROOT.TH1F('ctphiScintCal1d',
646  expRun + 'Scintillator BKLMHit1d phi-strip ctime distribution;' +
647  't(phi1D) - dt(layer) (ns)',
648  128, -0.5, 1023.5)
649 
650  self.hist_ctzScintCal1d = ROOT.TH1F('ctzScintCal1d',
651  expRun + 'Scintillator BKLMHit1d z-strip ctime distribution;' +
652  't(z1D) - dt(layer) (ns)',
653  128, -0.5, 1023.5)
654 
655  self.hist_ctScintCal1d = ROOT.TH1F('ctScintCal1d',
656  expRun + 'Scintillator BKLMHit1d x 2 calibrated average-time distribution;' +
657  '0.5*[t(phi1D) + t(z1D)] - dt(layer) (ns)',
658  128, -0.5, 1023.5)
659 
660  self.hist_dtScint1d = ROOT.TH1F('dtScint1d',
661  expRun + 'Scintillator BKLMHit1d x 2 time-difference distribution;' +
662  't(phi1D) - t(z1D) (ns)',
663  50, -100.0, 100.0)
664 
665  # Create the BKLMHit2d-related histograms
666 
667 
668  self.hist_nHit2d = ROOT.TH1F('NHit2d', expRun + '# of BKLMHit2ds', 50, -0.5, 49.5)
669 
670  self.hist_occupancyForwardXYPrompt = ROOT.TH2F('occupancyForwardXYPrompt',
671  expRun + 'Forward xy RPC occupancy for in-time hits;x(cm);y(cm)',
672  230, -345.0, 345.0, 230, -345.0, 345.0)
673 
674  self.hist_occupancyBackwardXYPrompt = ROOT.TH2F('occupancyBackwardXYPrompt',
675  expRun + 'Backward xy RPC occupancy for in-time hits;x(cm);y(cm)',
676  230, -345.0, 345.0, 230, -345.0, 345.0)
677 
678  self.hist_occupancyForwardXYBkgd = ROOT.TH2F('occupancyForwardXYBkgd',
679  expRun + 'Forward xy RPC occupancy for out-of-time hits;x(cm);y(cm)',
680  230, -345.0, 345.0, 230, -345.0, 345.0)
681 
682  self.hist_occupancyBackwardXYBkgd = ROOT.TH2F('occupancyBackwardXYBkgd',
683  expRun + 'Backward xy RPC occupancy for out-of-time hits;x(cm);y(cm)',
684  230, -345.0, 345.0, 230, -345.0, 345.0)
685 
686  self.hist_occupancyRZPrompt = ROOT.TH2F('occupancyRZPrompt',
687  expRun + 'layer-z occupancy for in-time hits;z(cm);layer',
688  48, -190.0, 290.0, 16, -0.5, 15.5)
689 
690  self.hist_occupancyZPrompt = ROOT.TH1F('occupancyZPrompt',
691  expRun + 'z occupancy for in-time hits;z(cm)',
692  48, -190.0, 290.0)
693 
694  self.hist_occupancyRPrompt = ROOT.TH1F('occupancyRPrompt',
695  expRun + 'layer occupancy for in-time hits;layer',
696  16, -0.5, 15.5)
697 
698  self.hist_occupancyRZBkgd = ROOT.TH2F('occupancyRZBkgd',
699  expRun + 'layer-z occupancy for out-of-time hits;z(cm);layer',
700  48, -190.0, 290.0, 16, -0.5, 15.5)
701 
702  self.hist_occupancyZBkgd = ROOT.TH1F('occupancyZBkgd',
703  expRun + 'z occupancy for out-of-time hits;z(cm)',
704  48, -190.0, 290.0)
705 
706  self.hist_occupancyRBkgd = ROOT.TH1F('occupancyRBkgd',
707  expRun + 'layer occupancy for out-of-time hits;layer',
708  16, -0.5, 15.5)
709 
710  self.hist_tRPCCal2d = ROOT.TH1F('tRPCCal2d',
711  expRun + 'RPC BKLMHit2d time distribution;' +
712  't(2D) - dt(layer) (ns)',
713  256, -0.5, 1023.5)
714 
715  self.hist_tRPCCal2dBySector = ROOT.TH2F('tRPCCal2dBySector',
716  expRun + 'RPC BKLMHit2d time distribution;' +
717  'sector # (0-7 = backward, 8-15 = forward);' +
718  't(2D) - dt(layer) (ns)',
719  16, -0.5, 15.5, 256, -0.5, 1023.5)
720 
721  self.hist_ctScintCal2d = ROOT.TH1F('ctScintCal2d',
722  expRun + 'Scint BKLMHit2d ctime distribution;' +
723  't(2D) - dt(layer) (ns)',
724  128, -0.5, 1023.5)
725 
726  self.hist_ctScintCal2dBySector = ROOT.TH2F('ctScintCal2dBySector',
727  expRun + 'Scint BKLMHit2d ctime distribution;' +
728  'sector # (0-7 = backward, 8-15 = forward);' +
729  't(2D) - dt(layer) (ns)',
730  16, -0.5, 15.5, 128, -0.5, 1023.5)
731 
732 
733  self.hist_tVsZFwd = ROOT.TProfile("tVsZForward",
734  expRun + 'RPC BKLMHit2d time vs z (forward);' +
735  'z (cm);' +
736  't(2D) - dt(layer) (ns)',
737  48, 0.0, 216.96)
738 
739  self.hist_tVsZBwd = ROOT.TProfile("tVsZBackward",
740  expRun + 'RPC BKLMHit2d time vs z (backward);' +
741  'z (cm);' +
742  't(2D) - dt(layer) (ns)',
743  48, 0.0, 216.96)
744  # Open the output PDF file for event displays
745 
746  if self.maxDisplays > 0:
747 
748  self.eventCanvas = ROOT.TCanvas("eventCanvas", self.eventPdfName, 3200, 1600)
749  title = '{0}['.format(self.eventPdfName)
750  self.eventCanvas.SaveAs(title)
751  self.eventCanvas.Clear()
752  self.eventCanvas.Divide(2, 1)
753 
754  # Create the boilerplate for the end- and side-views of the event display
755 
756 
757  self.bklmXY = []
758  r0 = 201.9 + 0.5 * 4.4 # cm
759  dr = 9.1 # cm
760  tan0 = math.tan(math.pi / 8.0)
761  g = ROOT.TGraph()
762  g.SetPoint(0, -200.0, 0.0)
763  g.SetPoint(1, +200.0, 0.0)
764  g.SetLineColor(19)
765  g.SetLineWidth(1)
766  g.SetLineStyle(3)
767  self.bklmXY.append(g)
768  g = ROOT.TGraph()
769  g.SetPoint(0, 0.0, -200.0)
770  g.SetPoint(1, 0.0, +200.0)
771  g.SetLineColor(19)
772  g.SetLineWidth(1)
773  g.SetLineStyle(3)
774  self.bklmXY.append(g)
775  g = ROOT.TGraph()
776  g.SetPoint(0, -5.0, 0.0)
777  g.SetPoint(1, +5.0, 0.0)
778  g.SetPoint(2, 0.0, 0.0)
779  g.SetPoint(3, 0.0, +5.0)
780  g.SetPoint(4, 0.0, -5.0)
781  g.SetLineColor(1)
782  g.SetLineWidth(1)
783  self.bklmXY.append(g)
784  for layer in range(0, 15):
785  r = r0 + layer * dr
786  x = r * tan0
787  g = ROOT.TGraph()
788  g.SetPoint(0, +r, -x)
789  g.SetPoint(1, +r, +x)
790  g.SetPoint(2, +x, +r)
791  g.SetPoint(3, -x, +r)
792  g.SetPoint(4, -r, +x)
793  g.SetPoint(5, -r, -x)
794  g.SetPoint(6, -x, -r)
795  g.SetPoint(7, +x, -r)
796  g.SetPoint(8, +r, -x)
797  if layer < 2:
798  g.SetLineColor(18)
799  else:
800  g.SetLineColor(17)
801  if (layer % 5) == 0:
802  g.SetLineStyle(3)
803  g.SetLineWidth(1)
804  self.bklmXY.append(g)
805 
806  self.bklmZY = []
807  rF = r0 + 14 * dr
808  x0 = r0 * tan0
809  z0 = 47.0 # cm
810  zL = 220.0 # cm
811  g = ROOT.TGraph()
812  g.SetPoint(0, -zL + z0 - 140.0, 0.0)
813  g.SetPoint(1, +zL + z0 + 70.0, 0.0)
814  g.SetLineColor(19)
815  g.SetLineWidth(1)
816  g.SetLineStyle(3)
817  self.bklmZY.append(g)
818  g = ROOT.TGraph()
819  g.SetPoint(0, 0.0, -315.0)
820  g.SetPoint(1, 0.0, +340.0)
821  g.SetLineColor(19)
822  g.SetLineWidth(1)
823  g.SetLineStyle(3)
824  self.bklmZY.append(g)
825  g = ROOT.TGraph()
826  g.SetPoint(0, -5.0, 0.0)
827  g.SetPoint(1, +5.0, 0.0)
828  g.SetPoint(2, 0.0, 0.0)
829  g.SetPoint(3, 0.0, +5.0)
830  g.SetPoint(4, 0.0, -5.0)
831  g.SetLineColor(1)
832  g.SetLineWidth(1)
833  self.bklmZY.append(g)
834  g = ROOT.TGraph()
835  g.SetPoint(0, -zL + z0, +x0)
836  g.SetPoint(1, -zL + z0, +r0)
837  g.SetLineColor(18)
838  g.SetLineWidth(1)
839  g.SetLineStyle(3)
840  self.bklmZY.append(g)
841  g = ROOT.TGraph()
842  g.SetPoint(0, -zL + z0, -x0)
843  g.SetPoint(1, -zL + z0, -r0)
844  g.SetLineColor(18)
845  g.SetLineWidth(1)
846  g.SetLineStyle(3)
847  self.bklmZY.append(g)
848  g = ROOT.TGraph()
849  g.SetPoint(0, +zL + z0, +x0)
850  g.SetPoint(1, +zL + z0, +r0)
851  g.SetLineColor(18)
852  g.SetLineWidth(1)
853  g.SetLineStyle(3)
854  self.bklmZY.append(g)
855  g = ROOT.TGraph()
856  g.SetPoint(0, +zL + z0, -x0)
857  g.SetPoint(1, +zL + z0, -r0)
858  g.SetLineColor(18)
859  g.SetLineWidth(1)
860  g.SetLineStyle(3)
861  self.bklmZY.append(g)
862  g = ROOT.TGraph()
863  g.SetPoint(0, -zL + z0, r0)
864  g.SetPoint(1, +zL + z0, r0)
865  g.SetPoint(2, +zL + z0, rF)
866  g.SetPoint(3, -zL + z0, rF)
867  g.SetPoint(4, -zL + z0, r0)
868  g.SetLineColor(18)
869  g.SetLineWidth(1)
870  self.bklmZY.append(g)
871  g = ROOT.TGraph()
872  g.SetPoint(0, -zL + z0, -r0)
873  g.SetPoint(1, +zL + z0, -r0)
874  g.SetPoint(2, +zL + z0, -rF)
875  g.SetPoint(3, -zL + z0, -rF)
876  g.SetPoint(4, -zL + z0, -r0)
877  g.SetLineColor(18)
878  g.SetLineWidth(1)
879  self.bklmZY.append(g)
880  g = ROOT.TGraph()
881  g.SetPoint(0, -zL + z0, -x0)
882  g.SetPoint(1, +zL + z0, -x0)
883  g.SetPoint(2, +zL + z0, +x0)
884  g.SetPoint(3, -zL + z0, +x0)
885  g.SetPoint(4, -zL + z0, -x0)
886  g.SetLineColor(18)
887  g.SetLineWidth(1)
888  self.bklmZY.append(g)
889 
890  def terminate(self):
891  """Handle job termination: draw histograms, close output files"""
892 
893  if self.maxDisplays > 0:
894  pdfNameLast = '{0}]'.format(self.eventPdfName)
895  self.eventCanvas.Print(pdfNameLast, self.lastTitle)
896 
897  for sectorFB in range(0, 16):
898  mappedScintSectorOccupancy = self.hist_mappedScintSectorOccupancy.GetBinContent(sectorFB + 1)
899  if mappedScintSectorOccupancy > 0:
900  for laneAxis in range(0, 42):
901  numerator = self.hist_mappedScintLaneAxisOccupancy.GetBinContent(sectorFB + 1, laneAxis + 1)
902  self.hist_mappedScintLaneAxisOccupancy.SetBinContent(
903  sectorFB + 1, laneAxis + 1, 100.0 * numerator / mappedScintSectorOccupancy)
904  mappedRPCSectorOccupancy = self.hist_mappedRPCSectorOccupancy.GetBinContent(sectorFB + 1)
905  if mappedRPCSectorOccupancy > 0:
906  for laneAxis in range(0, 42):
907  numerator = self.hist_mappedRPCLaneAxisOccupancy.GetBinContent(sectorFB + 1, laneAxis + 1)
908  self.hist_mappedRPCLaneAxisOccupancy.SetBinContent(
909  sectorFB + 1, laneAxis + 1, 100.0 * numerator / mappedRPCSectorOccupancy)
910  unmappedScintSectorOccupancy = self.hist_unmappedScintSectorOccupancy.GetBinContent(sectorFB + 1)
911  if unmappedScintSectorOccupancy > 0:
912  for laneAxis in range(0, 42):
913  numerator = self.hist_unmappedScintLaneAxisOccupancy.GetBinContent(sectorFB + 1, laneAxis + 1)
914  self.hist_unmappedScintLaneAxisOccupancy.SetBinContent(
915  sectorFB + 1, laneAxis + 1, 100.0 * numerator / unmappedScintSectorOccupancy)
916  unmappedRPCSectorOccupancy = self.hist_unmappedRPCSectorOccupancy.GetBinContent(sectorFB + 1)
917  if unmappedRPCSectorOccupancy > 0:
918  for laneAxis in range(0, 42):
919  numerator = self.hist_unmappedRPCLaneAxisOccupancy.GetBinContent(sectorFB + 1, laneAxis + 1)
920  self.hist_unmappedRPCLaneAxisOccupancy.SetBinContent(
921  sectorFB + 1, laneAxis + 1, 100.0 * numerator / unmappedRPCSectorOccupancy)
922  canvas = ROOT.TCanvas("canvas", self.pdfName, 1600, 1600)
923  title = '{0}['.format(self.pdfName)
924  canvas.SaveAs(title)
925  canvas.Clear()
926  canvas.GetPad(0).SetGrid(1, 1)
927  canvas.GetPad(0).Update()
928  self.hist_nDigit.Draw()
929  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nDigit.GetName()))
930  if self.verbosity > 0:
931  self.hist_nRawKLM.Draw()
932  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nRawKLM.GetName()))
933  self.hist_rawKLMnumEvents.Draw()
934  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMnumEvents.GetName()))
935  self.hist_rawKLMnumNodes.Draw()
936  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMnumNodes.GetName()))
937  self.hist_rawKLMnodeID.Draw("box")
938  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMnodeID.GetName()))
939  self.hist_rawKLMlaneFlag.Draw("box")
940  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMlaneFlag.GetName()))
941  # self.hist_rawKLMtdcExtraRPC.Draw("box")
942  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMtdcExtraRPC.GetName()))
943  # self.hist_rawKLMadcExtraRPC.Draw("box")
944  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMadcExtraRPC.GetName()))
945  # self.hist_rawKLMtdcExtraScint.Draw("box")
946  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMtdcExtraScint.GetName()))
947  # self.hist_rawKLMadcExtraScint.Draw("box")
948  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMadcExtraScint.GetName()))
949  if self.verbosity > 0:
950  self.hist_rawKLMsizeMultihit.Draw()
951  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMsizeMultihit.GetName()))
952  for dc in range(0, 16):
953  self.hist_rawKLMsizeByDCMultihit[dc].Draw()
954  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMsizeByDCMultihit[dc].GetName()))
955  self.hist_rawKLMsize.Draw()
956  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMsize.GetName()))
957  if self.verbosity > 0:
958  for dc in range(0, 16):
959  self.hist_rawKLMsizeByDC[dc].Draw()
960  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMsizeByDC[dc].GetName()))
961  for dc in range(0, 16):
962  self.hist_rawKLMchannelMultiplicity[dc].Draw("box")
963  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_rawKLMchannelMultiplicity[dc].GetName()))
964  # self.hist_mappedSectorOccupancy.Draw()
965  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedSectorOccupancy.GetName()))
966  # self.hist_unmappedSectorOccupancy.Draw()
967  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedSectorOccupancy.GetName()))
968  if self.verbosity > 0:
970  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCSectorOccupancy.GetName()))
972  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedRPCSectorOccupancy.GetName()))
974  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedScintSectorOccupancy.GetName()))
976  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintSectorOccupancy.GetName()))
977  canvas.Clear()
978  canvas.Divide(2, 1)
979  canvas.GetPad(0).SetGrid(1, 1)
980  canvas.GetPad(1).SetGrid(1, 1)
981  canvas.GetPad(2).SetGrid(1, 1)
982  # border of mapped z-readout channels for RPCs (axis=0)
983  borderRPC0x = [7.5, 20.5, 20.5, 7.5, 7.5]
984  borderRPC0y = [0.5, 0.5, 48.5, 48.5, 0.5]
985  borderRPC0yChimney = [0.5, 0.5, 34.5, 34.5, 0.5]
986  # border of mapped phi-readout channels for scints (axis=0)
987  borderScint0x = [0.5, 1.5, 1.5, 2.5, 2.5, 1.5, 1.5, 0.5, 0.5]
988  borderScint0y = [4.5, 4.5, 2.5, 2.5, 44.5, 44.5, 41.5, 41.5, 4.5]
989  # border of mapped phi-readout channels for RPCs (axis=1)
990  borderRPC1x = [7.5, 20.5, 20.5, 11.5, 11.5, 7.5, 7.5]
991  borderRPC1y = [0.5, 0.5, 48.5, 48.5, 36.5, 36.5, 0.5]
992  # border of mapped z-readout channels for scints (axis=1)
993  borderScint1x = [0.5, 2.5, 2.5, 0.5, 0.5]
994  borderScint1ay = [0.5, 0.5, 9.5, 9.5, 0.5]
995  borderScint1by = [15.5, 15.5, 60.5, 60.5, 15.5]
996  borderScint1xChimney = [0.5, 1.5, 1.5, 2.5, 2.5, 1.5, 1.5, 0.5, 0.5]
997  borderScint1ayChimney = [0.5, 0.5, 0.5, 0.5, 9.5, 9.5, 8.5, 8.5, 0.5]
998  borderScint1byChimney = [15.5, 15.5, 16.5, 16.5, 45.5, 45.5, 45.5, 45.5, 15.5]
999  # graphs of z-readout-channel borders for RPCs (axis=0)
1000  graphRPC0 = self.makeGraph(borderRPC0x, borderRPC0y)
1001  graphRPC0Chimney = self.makeGraph(borderRPC0x, borderRPC0yChimney)
1002  # graph of phi-readout-channel border for scints (axis=0)
1003  graphScint0 = self.makeGraph(borderScint0x, borderScint0y)
1004  # graph of phi-readout-channel border for RPCs (axis=1)
1005  graphRPC1 = self.makeGraph(borderRPC1x, borderRPC1y)
1006  # graphs of z-readout-channel borders for scints (axis=1)
1007  graphScint1a = self.makeGraph(borderScint1x, borderScint1ay)
1008  graphScint1b = self.makeGraph(borderScint1x, borderScint1by)
1009  graphScint1aChimney = self.makeGraph(borderScint1xChimney, borderScint1ayChimney)
1010  graphScint1bChimney = self.makeGraph(borderScint1xChimney, borderScint1byChimney)
1011  # labels for the above borders
1012  textRPC0 = self.makeText(6.8, 25.0, "RPC z")
1013  textScint0 = self.makeText(3.2, 25.0, "Scint #phi")
1014  textRPC1 = self.makeText(6.8, 25.0, "RPC #phi")
1015  textScint1 = self.makeText(3.2, 25.0, "Scint z")
1016  for sectorFB in range(0, 16):
1017  zmax = 1
1018  for lane in range(8, 21):
1019  for channel in range(0, 64):
1020  z = self.hist_mappedChannelOccupancyPrompt[sectorFB][0].GetBinContent(lane + 1, channel + 1)
1021  zmax = z if z > zmax else zmax
1022  self.hist_mappedChannelOccupancyPrompt[sectorFB][0].SetMaximum(zmax)
1023  canvas.cd(1)
1024  self.hist_mappedChannelOccupancyPrompt[sectorFB][0].Draw("colz") # z hits
1025  if sectorFB == 2:
1026  graphRPC0Chimney.Draw("L")
1027  graphScint1aChimney.Draw("L")
1028  graphScint1bChimney.Draw("L")
1029  else:
1030  graphRPC0.Draw("L")
1031  graphScint1a.Draw("L")
1032  graphScint1b.Draw("L")
1033  textRPC0.Draw()
1034  textScint1.Draw()
1035  zmax = 1
1036  for lane in range(8, 21):
1037  for channel in range(0, 64):
1038  z = self.hist_mappedChannelOccupancyPrompt[sectorFB][1].GetBinContent(lane + 1, channel + 1)
1039  zmax = z if z > zmax else zmax
1040  self.hist_mappedChannelOccupancyPrompt[sectorFB][1].SetMaximum(zmax)
1041  canvas.cd(2)
1042  self.hist_mappedChannelOccupancyPrompt[sectorFB][1].Draw("colz") # phi hits
1043  graphRPC1.Draw("L")
1044  graphScint0.Draw("L")
1045  textRPC1.Draw()
1046  textScint0.Draw()
1047  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedChannelOccupancyPrompt[sectorFB][0].GetName()))
1048  if self.verbosity > 0:
1049  for sectorFB in range(0, 16):
1050  zmax = 1
1051  for lane in range(8, 21):
1052  for channel in range(0, 64):
1053  z = self.hist_mappedChannelOccupancyBkgd[sectorFB][0].GetBinContent(lane + 1, channel + 1)
1054  zmax = z if z > zmax else zmax
1055  self.hist_mappedChannelOccupancyBkgd[sectorFB][0].SetMaximum(zmax)
1056  canvas.cd(1)
1057  self.hist_mappedChannelOccupancyBkgd[sectorFB][0].Draw("colz") # z hits
1058  if sectorFB == 2:
1059  graphRPC0Chimney.Draw("L")
1060  graphScint1aChimney.Draw("L")
1061  graphScint1bChimney.Draw("L")
1062  else:
1063  graphRPC0.Draw("L")
1064  graphScint1a.Draw("L")
1065  graphScint1b.Draw("L")
1066  textRPC0.Draw()
1067  textScint1.Draw()
1068  zmax = 1
1069  for lane in range(8, 21):
1070  for channel in range(0, 64):
1071  z = self.hist_mappedChannelOccupancyBkgd[sectorFB][1].GetBinContent(lane + 1, channel + 1)
1072  zmax = z if z > zmax else zmax
1073  self.hist_mappedChannelOccupancyBkgd[sectorFB][1].SetMaximum(zmax)
1074  canvas.cd(2)
1075  self.hist_mappedChannelOccupancyBkgd[sectorFB][1].Draw("colz") # phi hits
1076  graphRPC1.Draw("L")
1077  graphScint0.Draw("L")
1078  textRPC1.Draw()
1079  textScint0.Draw()
1080  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedChannelOccupancyBkgd[sectorFB][0].GetName()))
1081  for sectorFB in range(0, 16):
1082  canvas.cd(1)
1083  self.hist_unmappedChannelOccupancy[sectorFB][0].Draw("colz") # z hits
1084  if sectorFB == 2:
1085  graphRPC0Chimney.Draw("L")
1086  graphScint1aChimney.Draw("L")
1087  graphScint1bChimney.Draw("L")
1088  else:
1089  graphRPC0.Draw("L")
1090  graphScint1a.Draw("L")
1091  graphScint1b.Draw("L")
1092  textRPC0.Draw()
1093  textScint1.Draw()
1094  canvas.cd(2)
1095  self.hist_unmappedChannelOccupancy[sectorFB][1].Draw("colz") # phi hits
1096  graphRPC1.Draw("L")
1097  graphScint0.Draw("L")
1098  textRPC1.Draw()
1099  textScint0.Draw()
1100  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedChannelOccupancy[sectorFB][0].GetName()))
1101  canvas.Clear()
1102  canvas.Divide(1, 1)
1103  # self.hist_RPCTimeLowBitsBySector.Draw("box")
1104  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_RPCTimeLowBitsBySector.GetName()))
1105  timeFit = ROOT.TF1("timeFit", "gausn(0)+pol0(3)", 100.0, 500.0)
1106  timeFit.SetParName(0, "Area")
1107  timeFit.SetParName(1, "Mean")
1108  timeFit.SetParName(2, "Sigma")
1109  timeFit.SetParName(3, "Bkgd")
1110  timeFit.SetParLimits(2, 10.0, 80.0)
1111  n = self.hist_mappedRPCTime.GetEntries()
1112  timeFit.SetParameter(0, n)
1113  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1114  timeFit.SetParameter(1, 320.0)
1115  timeFit.SetParameter(2, 20.0)
1116  timeFit.SetParameter(3, 100.0)
1117  self.hist_mappedRPCTime.Fit("timeFit", "QR")
1118  self.hist_mappedRPCTime.Draw()
1119  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCTime.GetName()))
1120  self.hist_mappedRPCTimeBySector.Draw("box")
1121  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCTimeBySector.GetName()))
1122  n = self.hist_mappedRPCTimeCal.GetEntries()
1123  timeFit.SetParameter(0, n)
1124  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1125  timeFit.SetParameter(1, 320.0)
1126  timeFit.SetParameter(2, 20.0)
1127  timeFit.SetParameter(3, 100.0)
1128  self.hist_mappedRPCTimeCal.Fit("timeFit", "QR")
1129  self.hist_mappedRPCTimeCal.Draw()
1130  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCTimeCal.GetName()))
1131  self.hist_mappedRPCTimeCalBySector.Draw("box")
1132  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCTimeCalBySector.GetName()))
1133  if self.verbosity > 0:
1134  self.hist_unmappedRPCTime.Draw()
1135  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedRPCTime.GetName()))
1136  self.hist_unmappedRPCTimeBySector.Draw("box")
1137  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedRPCTimeBySector.GetName()))
1138  for sectorFB in range(0, 16):
1139  for layer in range(2, 15):
1140  n = self.hist_mappedRPCPhiTimePerLayer[sectorFB][layer].GetEntries()
1141  timeFit.SetParameter(0, n)
1142  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1143  timeFit.SetParameter(1, 320.0)
1144  timeFit.SetParameter(2, 20.0)
1145  timeFit.SetParameter(3, 100.0)
1146  self.hist_mappedRPCPhiTimePerLayer[sectorFB][layer].Fit("timeFit", "QR")
1147  self.hist_mappedRPCPhiTimePerLayer[sectorFB][layer].Draw()
1148  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCPhiTimePerLayer[sectorFB][layer].GetName()))
1149  for sectorFB in range(0, 16):
1150  for layer in range(2, 15):
1151  n = self.hist_mappedRPCZTimePerLayer[sectorFB][layer].GetEntries()
1152  timeFit.SetParameter(0, n)
1153  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1154  timeFit.SetParameter(1, 320.0)
1155  timeFit.SetParameter(2, 20.0)
1156  timeFit.SetParameter(3, 100.0)
1157  self.hist_mappedRPCZTimePerLayer[sectorFB][layer].Fit("timeFit", "QR")
1158  self.hist_mappedRPCZTimePerLayer[sectorFB][layer].Draw()
1159  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCZTimePerLayer[sectorFB][layer].GetName()))
1160  # self.hist_ScintTimeLowBitsBySector.Draw("box")
1161  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ScintTimeLowBitsBySector.GetName()))
1162  # self.hist_mappedScintTDC.Draw()
1163  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintTDC.GetName()))
1164  # self.hist_mappedScintTDCBySector.Draw("box")
1165  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintTDCBySector.GetName()))
1166  # self.hist_mappedScintTime.Draw()
1167  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintTime.GetName()))
1168  # self.hist_mappedScintTimeBySector.Draw("box")
1169  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintTimeBySector.GetName()))
1170  # self.hist_unmappedScintTime.Draw()
1171  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedScintTime.GetName()))
1172  # self.hist_unmappedScintTimeBySector.Draw("box")
1173  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedScintTimeBySector.GetName()))
1174  ctimeFit = ROOT.TF1("ctimeFit", "gausn(0)+pol0(3)", 300.0, 700.0)
1175  ctimeFit.SetParName(0, "Area")
1176  ctimeFit.SetParName(1, "Mean")
1177  ctimeFit.SetParName(2, "Sigma")
1178  ctimeFit.SetParName(3, "Bkgd")
1179  ctimeFit.SetParLimits(2, 10.0, 80.0)
1180  n = 32 * self.hist_mappedScintCtime.GetEntries()
1181  ctimeFit.SetParameter(0, n)
1182  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1183  ctimeFit.SetParameter(1, 450.0)
1184  ctimeFit.SetParameter(2, 40.0)
1185  ctimeFit.SetParameter(3, 10.0)
1186  self.hist_mappedScintCtime.Fit("ctimeFit", "QR")
1187  self.hist_mappedScintCtime.Draw()
1188  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtime.GetName()))
1189  self.hist_mappedScintCtimeBySector.Draw("box")
1190  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtimeBySector.GetName()))
1191  n = 32 * self.hist_mappedScintCtimeCal.GetEntries()
1192  ctimeFit.SetParameter(0, n)
1193  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1194  ctimeFit.SetParameter(1, 450.0)
1195  ctimeFit.SetParameter(2, 40.0)
1196  ctimeFit.SetParameter(3, 10.0)
1197  self.hist_mappedScintCtimeCal.Fit("ctimeFit", "QR")
1198  self.hist_mappedScintCtimeCal.Draw()
1199  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtimeCal.GetName()))
1200  self.hist_mappedScintCtimeCalBySector.Draw("box")
1201  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtimeCalBySector.GetName()))
1202  if self.verbosity > 0:
1203  self.hist_unmappedScintCtime.Draw()
1204  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedScintCtime.GetName()))
1205  self.hist_unmappedScintCtimeBySector.Draw("box")
1206  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_unmappedScintCtimeBySector.GetName()))
1207  for sectorFB in range(0, 16):
1208  for layer in range(0, 2):
1209  n = 32 * self.hist_mappedScintPhiCtimePerLayer[sectorFB][layer].GetEntries()
1210  ctimeFit.SetParameter(0, n)
1211  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1212  ctimeFit.SetParameter(1, 440.0)
1213  ctimeFit.SetParameter(2, 40.0)
1214  ctimeFit.SetParameter(3, 10.0)
1215  self.hist_mappedScintPhiCtimePerLayer[sectorFB][layer].Fit("ctimeFit", "QR")
1216  self.hist_mappedScintPhiCtimePerLayer[sectorFB][layer].Draw()
1217  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintPhiCtimePerLayer[sectorFB][layer].GetName()))
1218  for sectorFB in range(0, 16):
1219  for layer in range(0, 2):
1220  n = 32 * self.hist_mappedScintZCtimePerLayer[sectorFB][layer].GetEntries()
1221  ctimeFit.SetParameter(0, n)
1222  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1223  ctimeFit.SetParameter(1, 440.0)
1224  ctimeFit.SetParameter(2, 40.0)
1225  ctimeFit.SetParameter(3, 10.0)
1226  self.hist_mappedScintZCtimePerLayer[sectorFB][layer].Fit("ctimeFit", "QR")
1227  self.hist_mappedScintZCtimePerLayer[sectorFB][layer].Draw()
1228  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintZCtimePerLayer[sectorFB][layer].GetName()))
1229  self.hist_mappedRPCCtimeRangeBySector.Draw("box")
1230  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedRPCCtimeRangeBySector.GetName()))
1231  canvas.SetLogy(0)
1232  self.hist_tdcRangeRPC.Draw()
1233  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeRPC.GetName()))
1234  canvas.SetLogy(1)
1235  self.hist_tdcRangeRPC.Draw("")
1236  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeRPC.GetName()))
1237  canvas.SetLogy(0)
1238  self.hist_ctimeRangeRPC.Draw()
1239  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctimeRangeRPC.GetName()))
1240  canvas.SetLogy(1)
1241  self.hist_ctimeRangeRPC.Draw()
1242  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctimeRangeRPC.GetName()))
1243  canvas.SetLogy(0)
1244  canvas.SetLogz(0)
1245  self.hist_tdcRangeVsCtimeRangeRPC.Draw("BOX")
1246  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeVsCtimeRangeRPC.GetName()))
1247  canvas.SetLogz(1)
1248  self.hist_tdcRangeVsCtimeRangeRPC.Draw("COLZ")
1249  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeVsCtimeRangeRPC.GetName()))
1250  canvas.SetLogz(0)
1251  self.hist_tdcRangeVsTimeRPC.Draw("BOX")
1252  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeVsTimeRPC.GetName()))
1253  canvas.SetLogz(1)
1254  self.hist_tdcRangeVsTimeRPC.Draw("COLZ")
1255  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tdcRangeVsTimeRPC.GetName()))
1256  canvas.SetLogz(0)
1257  self.hist_ctimeRangeVsTimeRPC.Draw("BOX")
1258  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctimeRangeVsTimeRPC.GetName()))
1259  canvas.SetLogz(1)
1260  self.hist_ctimeRangeVsTimeRPC.Draw("COLZ")
1261  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctimeRangeVsTimeRPC.GetName()))
1262  canvas.SetLogz(0)
1263  self.hist_mappedScintCtimeRange.Draw()
1264  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtimeRange.GetName()))
1265  self.hist_mappedScintCtimeRangeBySector.Draw("box")
1266  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_mappedScintCtimeRangeBySector.GetName()))
1267 
1268  self.hist_nHit1d.Draw()
1269  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1d.GetName()))
1270  self.hist_n1dPhiZ.Draw("box")
1271  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_n1dPhiZ.GetName()))
1272  # self.hist_nHit1dRPCPrompt.Draw()
1273  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1dRPCPrompt.GetName()))
1274  # self.hist_nHit1dRPCBkgd.Draw()
1275  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1dRPCBkgd.GetName()))
1276  # self.hist_nHit1dScint.Draw()
1277  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1dScint.GetName()))
1278  # self.hist_nHit1dPrompt.Draw()
1279  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1dPrompt.GetName()))
1280  # self.hist_nHit1dBkgd.Draw()
1281  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit1dBkgd.GetName()))
1282  # self.hist_multiplicityPhiBySector.Draw("box")
1283  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_multiplicityPhiBySector.GetName()))
1284  # self.hist_multiplicityZBySector.Draw("box")
1285  # canvas.Print(self.pdfName, "Title:{0}".format(self.hist_multiplicityZBySector.GetName()))
1286  n = self.hist_tphiRPCCal1d.GetEntries()
1287  timeFit.SetParameter(0, n)
1288  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1289  timeFit.SetParameter(1, 330.0)
1290  timeFit.SetParameter(2, 15.0)
1291  timeFit.SetParameter(3, 100.0)
1292  self.hist_tphiRPCCal1d.Fit("timeFit", "QR")
1293  self.hist_tphiRPCCal1d.Draw()
1294  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tphiRPCCal1d.GetName()))
1295  n = self.hist_tzRPCCal1d.GetEntries()
1296  timeFit.SetParameter(0, n)
1297  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1298  timeFit.SetParameter(1, 330.0)
1299  timeFit.SetParameter(2, 15.0)
1300  timeFit.SetParameter(3, 100.0)
1301  self.hist_tzRPCCal1d.Fit("timeFit", "QR")
1302  self.hist_tzRPCCal1d.Draw()
1303  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tzRPCCal1d.GetName()))
1304  if self.verbosity > 0:
1305  n = self.hist_tRPCCal1d.GetEntries()
1306  timeFit.SetParameter(0, n)
1307  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1308  timeFit.SetParameter(1, 330.0)
1309  timeFit.SetParameter(2, 15.0)
1310  timeFit.SetParameter(3, 100.0)
1311  self.hist_tRPCCal1d.Fit("timeFit", "QR")
1312  self.hist_tRPCCal1d.Draw()
1313  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tRPCCal1d.GetName()))
1314  self.hist_dtRPC1d.Draw()
1315  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_dtRPC1d.GetName()))
1316  n = 32 * self.hist_ctphiScintCal1d.GetEntries()
1317  ctimeFit.SetParameter(0, n)
1318  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1319  ctimeFit.SetParameter(1, 520.0)
1320  ctimeFit.SetParameter(2, 40.0)
1321  ctimeFit.SetParameter(3, 10.0)
1322  self.hist_ctphiScintCal1d.Fit("ctimeFit", "QR")
1323  self.hist_ctphiScintCal1d.Draw()
1324  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctphiScintCal1d.GetName()))
1325  n = 32 * self.hist_ctzScintCal1d.GetEntries()
1326  ctimeFit.SetParameter(0, n)
1327  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1328  ctimeFit.SetParameter(1, 520.0)
1329  ctimeFit.SetParameter(2, 40.0)
1330  ctimeFit.SetParameter(3, 10.0)
1331  self.hist_ctzScintCal1d.Fit("ctimeFit", "QR")
1332  self.hist_ctzScintCal1d.Draw()
1333  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctzScintCal1d.GetName()))
1334  if self.verbosity > 0:
1335  n = 32 * self.hist_ctScintCal1d.GetEntries()
1336  ctimeFit.SetParameter(0, n)
1337  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1338  ctimeFit.SetParameter(1, 520.0)
1339  ctimeFit.SetParameter(2, 40.0)
1340  ctimeFit.SetParameter(3, 10.0)
1341  self.hist_ctScintCal1d.Fit("ctimeFit", "QR")
1342  self.hist_ctScintCal1d.Draw()
1343  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctScintCal1d.GetName()))
1344  self.hist_dtScint1d.Draw()
1345  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_dtScint1d.GetName()))
1346 
1347  self.hist_nHit2d.Draw()
1348  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_nHit2d.GetName()))
1349  self.hist_occupancyForwardXYPrompt.Draw("colz")
1350  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyForwardXYPrompt.GetName()))
1351  self.hist_occupancyBackwardXYPrompt.Draw("colz")
1352  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyBackwardXYPrompt.GetName()))
1353  if self.verbosity > 0:
1354  self.hist_occupancyForwardXYBkgd.Draw("colz")
1355  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyForwardXYBkgd.GetName()))
1356  self.hist_occupancyBackwardXYBkgd.Draw("colz")
1357  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyBackwardXYBkgd.GetName()))
1358  self.hist_occupancyRZPrompt.Draw("colz")
1359  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyRZPrompt.GetName()))
1360  self.hist_occupancyZPrompt.Draw()
1361  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyZPrompt.GetName()))
1362  self.hist_occupancyRPrompt.Draw()
1363  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyRPrompt.GetName()))
1364  self.hist_occupancyRZBkgd.Draw("colz")
1365  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyRZBkgd.GetName()))
1366  self.hist_occupancyZBkgd.Draw()
1367  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyZBkgd.GetName()))
1368  self.hist_occupancyRBkgd.Draw()
1369  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_occupancyRBkgd.GetName()))
1370  self.hist_tVsZFwd.Draw()
1371  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tVsZFwd.GetName()))
1372  self.hist_tVsZBwd.Draw()
1373  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tVsZFwd.GetName()))
1374  timeFit.SetParameter(0, self.hist_tRPCCal2d.GetEntries())
1375  n = self.hist_tRPCCal2d.GetEntries()
1376  timeFit.SetParameter(0, n)
1377  timeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1378  timeFit.SetParameter(1, 320.0)
1379  timeFit.SetParameter(2, 20.0)
1380  timeFit.SetParameter(3, 100.0)
1381  self.hist_tRPCCal2d.Fit("timeFit", "QR")
1382  self.hist_tRPCCal2d.Draw()
1383  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tRPCCal2d.GetName()))
1384  self.hist_tRPCCal2dBySector.Draw("box")
1385  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_tRPCCal2dBySector.GetName()))
1386  n = 32 * self.hist_ctScintCal2d.GetEntries()
1387  ctimeFit.SetParameter(0, n)
1388  ctimeFit.SetParLimits(0, 0.2 * n, 5.0 * n)
1389  ctimeFit.SetParameter(1, 500.0)
1390  ctimeFit.SetParameter(2, 40.0)
1391  ctimeFit.SetParameter(3, 10.0)
1392  self.hist_ctScintCal2d.Fit("ctimeFit", "QR")
1393  self.hist_ctScintCal2d.Draw()
1394  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctScintCal2d.GetName()))
1395  self.hist_ctScintCal2dBySector.Draw("box")
1396  canvas.Print(self.pdfName, "Title:{0}".format(self.hist_ctScintCal2dBySector.GetName()))
1397  # the last page of the PDF file has to have "]" at the end of the PDF filename
1398  pdfNameLast = '{0}]'.format(self.pdfName)
1399  canvas.Print(pdfNameLast, "Title:{0}".format(self.hist_ctScintCal2dBySector.GetName()))
1400  self.histogramFile.Write()
1401  self.histogramFile.Close()
1402  print('Goodbye')
1403 
1404  def beginRun(self):
1405  """Handle begin of run: print diagnostic message"""
1406  EventMetaData = Belle2.PyStoreObj('EventMetaData')
1407  print('beginRun', EventMetaData.getRun())
1408 
1409  def endRun(self):
1410  """Handle end of run: print diagnostic message"""
1411  EventMetaData = Belle2.PyStoreObj('EventMetaData')
1412  print('endRun', EventMetaData.getRun())
1413 
1414  def event(self):
1415  """Process one event: fill histograms, (optionally) draw event display"""
1416 
1417  self.eventCounter += 1
1418  EventMetaData = Belle2.PyStoreObj('EventMetaData')
1419  event = EventMetaData.getEvent()
1420  rawklms = Belle2.PyStoreArray('RawKLMs')
1421  digits = Belle2.PyStoreArray('BKLMDigits')
1422  hit1ds = Belle2.PyStoreArray('BKLMHit1ds')
1423  hit2ds = Belle2.PyStoreArray('BKLMHit2ds')
1424  self.hist_nRawKLM.Fill(len(rawklms))
1425  self.hist_nDigit.Fill(len(digits))
1426  self.hist_nHit1d.Fill(len(hit1ds))
1427  self.hist_nHit2d.Fill(len(hit2ds))
1428 
1429  tCal1d = []
1430  for hit1d in hit1ds:
1431  tCal1d.append(hit1d.getTime())
1432  tCal2d = []
1433  for hit2d in hit2ds:
1434  tCal2d.append(hit2d.getTime())
1435 
1436  countAllMultihit = 0
1437  countAll = 0
1438  count = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
1439  rawFb = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1440  rawSector = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1441  rawLayer = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1442  rawPlane = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1443  rawStrip = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1444  rawCtime = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
1445  for copper in range(0, len(rawklms)):
1446  rawklm = rawklms[copper]
1447  self.hist_rawKLMnumEvents.Fill(rawklm.GetNumEvents())
1448  self.hist_rawKLMnumNodes.Fill(rawklm.GetNumNodes())
1449  if rawklm.GetNumEntries() != 1:
1450  print('##0 Event', event, 'copper', copper, ' getNumEntries=', rawklm.GetNumEntries())
1451  continue
1452  nodeID = rawklm.GetNodeID(0) - self.BKLM_ID
1453  if nodeID >= self.EKLM_ID - self.BKLM_ID:
1454  nodeID = nodeID - (self.EKLM_ID - self.BKLM_ID) + 4
1455  self.hist_rawKLMnodeID.Fill(nodeID, copper)
1456  if (nodeID < 0) or (nodeID > 4): # skip EKLM nodes
1457  continue
1458  trigCtime = (rawklm.GetTTCtime(0) & 0x7ffffff) << 3 # (ns)
1459  revo9time = trigCtime - 0x3b0
1460  for finesse in range(0, 4):
1461  dc = (finesse << 2) + copper
1462  nWords = rawklm.GetDetectorNwords(0, finesse)
1463  self.hist_rawKLMsizeByDCMultihit[dc].Fill(nWords)
1464  if nWords <= 0:
1465  continue
1466  countAllMultihit = countAllMultihit + nWords
1467  bufSlot = rawklm.GetDetectorBuffer(0, finesse)
1468  lastWord = bufSlot[nWords - 1]
1469  if lastWord & 0xffff != 0:
1470  print("##1 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'lastWord=', hex(lastWord))
1471  if (nWords % 2) == 0:
1472  print("##2 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'should be odd -- skipping')
1473  continue
1474  if int(self.exp) != 3: # revo9time was not stored in the last word of the data-packet list?
1475  revo9time = ((lastWord >> 16) << 3) & 0xffff
1476  dt = (trigCtime - revo9time) & 0x3ff
1477  if dt >= 0x200:
1478  dt -= 0x400
1479  self.hist_trigCtimeVsTrigRevo9time.Fill(dt)
1480  countAll += 1
1481  count[dc] += 1
1482  sectorFB = self.dcToSectorFB[dc]
1483  n = nWords >> 1 # number of Data-Concentrator data packets
1484  channelMultiplicity = {}
1485  minRPCCtime = 99999
1486  maxRPCCtime = 0
1487  minRPCtdc = 99999
1488  maxRPCtdc = 0
1489  minScintCtime = 99999
1490  maxScintCtime = 0
1491  # first pass over this DC: determine per-channel multiplicities, event time ranges, and
1492  # fill dictionaries for accessing RawKLM hit information from BLKMHit1ds and BKLMHit2ds
1493  for j in range(0, n):
1494  word0 = bufSlot[j * 2]
1495  word1 = bufSlot[j * 2 + 1]
1496  ctime = word0 & 0xffff
1497  channel = (word0 >> 16) & 0x7f
1498  axis = (word0 >> 23) & 0x01
1499  lane = (word0 >> 24) & 0x1f # 1..2 for scints, 8..20 for RPCs (=readout-board slot - 7)
1500  flag = (word0 >> 30) & 0x03 # identifies scintillator or RPC hit
1501  adc = word1 & 0x0fff
1502  tdc = (word1 >> 16) & 0x07ff
1503  isRPC = (flag == 1)
1504  isScint = (flag == 2)
1505  laneAxisChannel = (word0 >> 16) & 0x1fff
1506  if laneAxisChannel not in channelMultiplicity:
1507  countAll = countAll + 2
1508  count[dc] = count[dc] + 2
1509  channelMultiplicity[laneAxisChannel] = 0
1510  channelMultiplicity[laneAxisChannel] += 1
1511  if isRPC:
1512  if ctime < minRPCCtime:
1513  minRPCCtime = ctime
1514  if ctime > maxRPCCtime:
1515  maxRPCCtime = ctime
1516  if tdc < minRPCtdc:
1517  minRPCtdc = tdc
1518  if tdc > maxRPCtdc:
1519  maxRPCtdc = tdc
1520  elif isScint:
1521  if int(self.exp) <= 3: # fix the ctime for old SCROD firmware
1522  trigCtx = trigCtime >> 3
1523  ctime = trigCtx - ((trigCtx - ctime) << 2)
1524  if ctime < minScintCtime:
1525  minScintCtime = ctime
1526  if ctime > maxScintCtime:
1527  maxScintCtime = ctime
1528  electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
1529  if electId in self.electIdToModuleId:
1530  moduleId = self.electIdToModuleId[electId]
1531  fb = (moduleId & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
1532  sector = (moduleId & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
1533  layer = (moduleId & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
1534  plane = (moduleId & self.BKLM_PLANE_MASK) >> self.BKLM_PLANE_BIT
1535  strip = (moduleId & self.BKLM_STRIP_MASK) >> self.BKLM_STRIP_BIT
1536  rawFb[dc].append(fb)
1537  rawSector[dc].append(sector)
1538  rawLayer[dc].append(layer)
1539  rawPlane[dc].append(plane)
1540  rawStrip[dc].append(strip)
1541  rawCtime[dc].append(ctime)
1542  else:
1543  rawFb[dc].append(-1)
1544  rawSector[dc].append(-1)
1545  rawLayer[dc].append(-1)
1546  rawPlane[dc].append(-1)
1547  rawStrip[dc].append(-1)
1548  rawCtime[dc].append(-1)
1549  tdcRangeRPC = maxRPCtdc - minRPCtdc # (ns)
1550  ctimeRangeRPC = (maxRPCCtime - minRPCCtime) << 3 # (ns)
1551  ctimeRangeScint = (maxScintCtime - minScintCtime) << 3 # (ns)
1552  if n > 1:
1553  if maxRPCCtime > 0:
1554  self.hist_mappedRPCCtimeRangeBySector.Fill(sectorFB, ctimeRangeRPC)
1555  if maxScintCtime > 0:
1556  self.hist_mappedScintCtimeRange.Fill(ctimeRangeScint)
1557  self.hist_mappedScintCtimeRangeBySector.Fill(sectorFB, ctimeRangeScint)
1558  # second pass over this DC's hits: histogram everything
1559  for j in range(0, n):
1560  word0 = bufSlot[j * 2]
1561  word1 = bufSlot[j * 2 + 1]
1562  ctime = word0 & 0xffff
1563  channel = (word0 >> 16) & 0x7f
1564  axis = (word0 >> 23) & 0x01
1565  lane = (word0 >> 24) & 0x1f # 1..2 for scints, 8..20 for RPCs (=readout-board slot - 7)
1566  flag = (word0 >> 30) & 0x03 # 1 for RPCs, 2 for scints
1567  electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
1568  adc = word1 & 0x0fff
1569  tdc = (word1 >> 16) & 0x07ff
1570  tdcExtra = (word1 >> 27) & 0x1f
1571  adcExtra = (word1 >> 12) & 0x0f
1572  isRPC = (flag == 1)
1573  isScint = (flag == 2)
1574  laneAxis = axis if ((lane < 1) or (lane > 20)) else ((lane << 1) + axis)
1575  laneAxisChannel = (word0 >> 16) & 0x1fff
1576  multiplicity = channelMultiplicity[laneAxisChannel]
1577  if multiplicity > 1: # histogram only if 2+ entries in the same channel
1578  self.hist_rawKLMchannelMultiplicity[dc].Fill(multiplicity, laneAxis)
1579  self.hist_rawKLMchannelMultiplicityFine[dc].Fill(multiplicity, laneAxisChannel)
1580  if (self.singleEntry == 1 and multiplicity > 1) or (self.singleEntry == 2 and multiplicity == 1):
1581  continue
1582  self.hist_rawKLMlaneFlag.Fill(flag, lane)
1583  if isRPC:
1584  self.hist_rawKLMtdcExtraRPC.Fill(sectorFB, tdcExtra)
1585  self.hist_rawKLMadcExtraRPC.Fill(sectorFB, adcExtra)
1586  elif isScint:
1587  self.hist_rawKLMtdcExtraScint.Fill(sectorFB, tdcExtra)
1588  self.hist_rawKLMadcExtraScint.Fill(sectorFB, adcExtra)
1589  if int(self.exp) <= 3: # fix the ctime for old SCROD firmware
1590  trigCtx = trigCtime >> 3
1591  ctime = trigCtx - ((trigCtx - ctime) << 2)
1592  t = (tdc - trigCtime) & 0x03ff # in ns, range is 0..1023
1593  ct = (ctime << 3) - (trigCtime & 0x7fff8) # in ns, range is only 8 bits in SCROD (??)
1594  ct = ct & 0x3ff
1595  if electId in self.electIdToModuleId: # mapped-channel histograms
1596  self.hist_mappedSectorOccupancyMultihit.Fill(sectorFB)
1597  if channelMultiplicity[laneAxisChannel] == 1:
1598  self.hist_mappedSectorOccupancy.Fill(sectorFB)
1599  if isRPC:
1600  self.hist_RPCTimeLowBitsBySector.Fill(sectorFB, (tdc & 3))
1601  tCal = int(tdc - trigCtime - self.t0RPC[axis][sectorFB][lane-6]) & 0x03ff # in ns, range 0..1023
1602  if j == 0:
1603  self.hist_tdcRangeRPC.Fill(tdcRangeRPC)
1604  self.hist_ctimeRangeRPC.Fill(ctimeRangeRPC)
1605  self.hist_tdcRangeVsCtimeRangeRPC.Fill(tdcRangeRPC, ctimeRangeRPC)
1606  self.hist_tdcRangeVsTimeRPC.Fill(tCal, tdcRangeRPC)
1607  self.hist_ctimeRangeVsTimeRPC.Fill(tCal, ctimeRangeRPC)
1608  if abs(tCal - self.t0Cal) < 50:
1609  self.hist_mappedChannelOccupancyPrompt[sectorFB][axis].Fill(lane, channel)
1610  else:
1611  self.hist_mappedChannelOccupancyBkgd[sectorFB][axis].Fill(lane, channel)
1612  self.hist_mappedRPCSectorOccupancy.Fill(sectorFB)
1613  self.hist_mappedRPCLaneAxisOccupancy.Fill(sectorFB, laneAxis)
1614  self.hist_mappedRPCTime.Fill(t)
1615  self.hist_mappedRPCTimeCal.Fill(tCal)
1616  self.hist_mappedRPCTimeBySector.Fill(sectorFB, t)
1617  self.hist_mappedRPCTimeCalBySector.Fill(sectorFB, tCal)
1618  if axis == 0:
1619  self.hist_mappedRPCZTimePerLayer[sectorFB][lane - 6].Fill(t)
1620  else:
1621  self.hist_mappedRPCPhiTimePerLayer[sectorFB][lane - 6].Fill(t)
1622  elif isScint:
1623  self.hist_ScintTimeLowBitsBySector.Fill(sectorFB, (tdc & 3))
1624  ctCal = int((ctime << 3) - trigCtime - self.ct0Scint[1-axis][sectorFB][lane-1]) & 0x03ff # in ns
1625  if abs(ctCal - self.ct0Cal) < 50:
1626  self.hist_mappedChannelOccupancyPrompt[sectorFB][1 - axis].Fill(lane, channel)
1627  else:
1628  self.hist_mappedChannelOccupancyBkgd[sectorFB][1 - axis].Fill(lane, channel)
1629  self.hist_mappedScintSectorOccupancy.Fill(sectorFB)
1630  self.hist_mappedScintLaneAxisOccupancy.Fill(sectorFB, laneAxis)
1631  self.hist_mappedScintTime.Fill(t & 0x1f)
1632  self.hist_mappedScintTimeBySector.Fill(sectorFB, t & 0x1f)
1633  self.hist_mappedScintTDC.Fill(tdc)
1634  self.hist_mappedScintTDCBySector.Fill(sectorFB, tdc)
1635  self.hist_mappedScintCtime.Fill(ct)
1636  self.hist_mappedScintCtimeBySector.Fill(sectorFB, ct)
1637  self.hist_mappedScintCtimeCal.Fill(ctCal)
1638  self.hist_mappedScintCtimeCalBySector.Fill(sectorFB, ctCal)
1639  if axis == 1:
1640  self.hist_mappedScintZCtimePerLayer[sectorFB][lane - 1].Fill(ct)
1641  else:
1642  self.hist_mappedScintPhiCtimePerLayer[sectorFB][lane - 1].Fill(ct)
1643  else: # unmapped-channel histograms
1644  self.hist_unmappedSectorOccupancyMultihit.Fill(sectorFB)
1645  if channelMultiplicity[laneAxisChannel] == 1:
1646  self.hist_unmappedSectorOccupancy.Fill(sectorFB)
1647  if isRPC:
1648  self.hist_unmappedChannelOccupancy[sectorFB][axis].Fill(lane, channel)
1649  self.hist_RPCTimeLowBitsBySector.Fill(sectorFB, (tdc & 3))
1650  self.hist_unmappedRPCSectorOccupancy.Fill(sectorFB)
1651  self.hist_unmappedRPCLaneAxisOccupancy.Fill(sectorFB, laneAxis)
1652  self.hist_unmappedRPCTime.Fill(t)
1653  self.hist_unmappedRPCTimeBySector.Fill(sectorFB, t)
1654  elif isScint:
1655  self.hist_unmappedChannelOccupancy[sectorFB][1 - axis].Fill(lane, channel)
1656  self.hist_ScintTimeLowBitsBySector.Fill(sectorFB, (tdc & 3))
1657  self.hist_unmappedScintSectorOccupancy.Fill(sectorFB)
1658  self.hist_unmappedScintLaneAxisOccupancy.Fill(sectorFB, laneAxis)
1659  self.hist_unmappedScintTime.Fill(t & 0x1f)
1660  self.hist_unmappedScintTimeBySector.Fill(sectorFB, t & 0x1f)
1661  self.hist_unmappedScintCtime.Fill(ct)
1662  self.hist_unmappedScintCtimeBySector.Fill(sectorFB, ct)
1663  self.hist_rawKLMsizeByDC[dc].Fill(count[dc])
1664  self.hist_rawKLMsizeMultihit.Fill(countAllMultihit)
1665  self.hist_rawKLMsize.Fill(countAll)
1666 
1667  # Process the BKLMHit1ds
1668 
1669  cosine = [0, 0, 0, 0, 0, 0, 0, 0]
1670  sine = [0, 0, 0, 0, 0, 0, 0, 0]
1671  for sector in range(0, 8):
1672  phi = math.pi * sector / 4
1673  cosine[sector] = math.cos(phi)
1674  sine[sector] = math.sin(phi)
1675  zyList = [[], [], [], [], [], [], [], []]
1676  xyList = [[], [], [], [], [], [], [], []]
1677  r0 = 201.9 + 0.5 * 4.4 # cm
1678  dr = 9.1 # cm
1679  z0 = 47.0 # cm
1680  dzScint = 4.0 # cm
1681  dzRPC = 4.52 # cm
1682  nPhiStrips = [37, 42, 36, 36, 36, 36, 48, 48, 48, 48, 48, 48, 48, 48, 48]
1683  dPhiStrips = [4.0, 4.0, 4.9, 5.11, 5.32, 5.53, 4.3, 4.46, 4.62, 4.77, 4.93, 5.09, 5.25, 5.4, 5.56]
1684  scintFlip = [[[-1, 1], [-1, 1], [1, 1], [1, -1], [1, -1], [1, -1], [-1, -1], [-1, 1]],
1685  [[1, -1], [1, -1], [1, 1], [-1, 1], [-1, 1], [-1, 1], [-1, -1], [1, -1]]]
1686  promptColor = 3
1687  bkgdColor = 2
1688  phiTimes = {}
1689  zTimes = [{}, {}, {}, {}, {}, {}, {}, {}, {}, {}, {}, {}, {}, {}, {}]
1690  nphihits = 0
1691  nzhits = 0
1692  nRPCPrompt = 0
1693  nRPCBkgd = 0
1694  nScint = 0
1695  for hit1d in hit1ds:
1696  key = hit1d.getModuleID()
1697  fb = (key & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
1698  sector = (key & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
1699  layer = (key & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
1700  plane = (key & self.BKLM_PLANE_MASK) >> self.BKLM_PLANE_BIT
1701  stripMin = (key & self.BKLM_STRIP_MASK) >> self.BKLM_STRIP_BIT
1702  stripMax = (key & self.BKLM_MAXSTRIP_MASK) >> self.BKLM_MAXSTRIP_BIT
1703  sectorFB = sector if fb == 0 else sector + 8
1704  if self.legacyTimes:
1705  dc = self.sectorFBToDC[sectorFB]
1706  copper = dc & 0x03
1707  finesse = dc >> 2
1708  n = rawklms[copper].GetDetectorNwords(0, finesse) >> 1
1709  trigCtime = (rawklms[copper].GetTTCtime(0) & 0x07ffffff) << 3
1710  tCal = -1
1711  ctDiffMax = 99999
1712  for j in range(0, n):
1713  if layer != rawLayer[dc][j]:
1714  continue
1715  if sector != rawSector[dc][j]:
1716  continue
1717  if fb != rawFb[dc][j]:
1718  continue
1719  if plane != rawPlane[dc][j]:
1720  continue
1721  strip = rawStrip[dc][j]
1722  if strip < stripMin:
1723  continue
1724  if strip > stripMax:
1725  continue
1726  if layer < 2: # it's a scint layer
1727  ctime = rawCtime[dc][j] << 3
1728  ct = ctime - trigCtime - self.ct0Scint[plane][sectorFB][layer]
1729  ctTrunc = int(ct) & 0x3ff
1730  if abs(ctTrunc - self.ct0Cal) < ctDiffMax:
1731  ctDiffMax = int(abs(ctTrunc - self.ct0Cal))
1732  tCal = ct
1733  if ctDiffMax == 0:
1734  break
1735  else: # it's an RPC layer
1736  tCal = tCal = hit1d.getTime() - trigCtime - self.t0RPC[plane][sectorFB][layer]
1737  break
1738  else:
1739  if layer < 2:
1740  tCal = hit1d.getTime() - self.ct0Scint[plane][sectorFB][layer]
1741  else:
1742  tCal = hit1d.getTime() - self.t0RPC[plane][sectorFB][layer]
1743  tCalTrunc = int(tCal) & 0x3ff
1744 
1745  if self.view == 1:
1746  r = r0 + layer * dr
1747  yA = r
1748  zA = 500
1749  xB = 500
1750  yB = 500
1751  stripAverage = (stripMin + stripMax) * 0.5
1752  isPrompt = False
1753  if layer < 2:
1754  nScint += 1
1755  isPrompt = (abs(tCalTrunc - self.ct0Cal1d) < 50)
1756  if plane == 0:
1757  if fb == 0:
1758  zA = z0 - stripAverage * dzScint
1759  else:
1760  zA = z0 + stripAverage * dzScint
1761  else:
1762  h = ((stripAverage - 0.5 * nPhiStrips[layer]) * dPhiStrips[layer]) * scintFlip[fb][sector][layer]
1763  xB = r * cosine[sector] - h * sine[sector]
1764  yB = r * sine[sector] + h * cosine[sector]
1765  else:
1766  isPrompt = (abs(tCalTrunc - self.t0Cal1d) < 50)
1767  if plane == 0:
1768  if fb == 0:
1769  zA = z0 - stripAverage * dzRPC
1770  else:
1771  zA = z0 + stripAverage * dzRPC
1772  else:
1773  h = ((stripAverage - 0.5 * nPhiStrips[layer]) * dPhiStrips[layer]) # * rpcFlip[fb][sector]
1774  xB = r * cosine[sector] - h * sine[sector]
1775  yB = r * sine[sector] + h * cosine[sector]
1776  if abs(tCalTrunc - self.t0Cal) < 50:
1777  nRPCPrompt += 1
1778  if plane == 1:
1779  self.hist_multiplicityPhiBySector.Fill(sectorFB, stripMax - stripMin + 1)
1780  else:
1781  self.hist_multiplicityZBySector.Fill(sectorFB, stripMax - stripMin + 1)
1782  else:
1783  nRPCBkgd += 1
1784  if plane == 1:
1785  nphihits += 1
1786  phiTimes[key] = tCal
1787  if layer < 2:
1788  self.hist_ctphiScintCal1d.Fill(tCalTrunc)
1789  else:
1790  self.hist_tphiRPCCal1d.Fill(tCalTrunc)
1791  else:
1792  nzhits += 1
1793  zTimes[layer][key] = tCal
1794  if layer < 2:
1795  self.hist_ctzScintCal1d.Fill(tCalTrunc)
1796  else:
1797  self.hist_tzRPCCal1d.Fill(tCalTrunc)
1798  # Add the hit to the event-display TGraph list (perhaps)
1799  if (self.view == 1) and (self.eventDisplays < self.maxDisplays):
1800  if zA != 500:
1801  gZY = ROOT.TGraph()
1802  gZY.SetPoint(0, zA - 1.0, yA - 1.0)
1803  gZY.SetPoint(1, zA - 1.0, yA + 1.0)
1804  gZY.SetPoint(2, zA + 1.0, yA + 1.0)
1805  gZY.SetPoint(3, zA + 1.0, yA - 1.0)
1806  gZY.SetPoint(4, zA - 1.0, yA - 1.0)
1807  gZY.SetLineWidth(1)
1808  if isPrompt:
1809  gZY.SetLineColor(promptColor)
1810  else:
1811  gZY.SetLineColor(bkgdColor)
1812  zyList[sector].append(gZY)
1813  if xB != 500:
1814  gXY = ROOT.TGraph()
1815  gXY.SetPoint(0, xB - 1.0, yB - 1.0)
1816  gXY.SetPoint(1, xB - 1.0, yB + 1.0)
1817  gXY.SetPoint(2, xB + 1.0, yB + 1.0)
1818  gXY.SetPoint(3, xB + 1.0, yB - 1.0)
1819  gXY.SetPoint(4, xB - 1.0, yB - 1.0)
1820  gXY.SetLineWidth(1)
1821  if isPrompt:
1822  gXY.SetLineColor(promptColor)
1823  else:
1824  gXY.SetLineColor(bkgdColor)
1825  xyList[sector].append(gXY)
1826  self.hist_nHit1dRPCPrompt.Fill(nRPCPrompt)
1827  self.hist_nHit1dRPCBkgd.Fill(nRPCBkgd)
1828  self.hist_nHit1dScint.Fill(nScint)
1829  if nRPCPrompt > 2:
1830  self.hist_nHit1dPrompt.Fill(nScint + nRPCBkgd + nRPCPrompt)
1831  else:
1832  self.hist_nHit1dBkgd.Fill(nScint + nRPCBkgd + nRPCPrompt)
1833  self.hist_n1dPhiZ.Fill(nphihits, nzhits)
1834  for phiKey in phiTimes:
1835  mphi = phiKey & self.BKLM_MODULEID_MASK
1836  layer = (mphi & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
1837  sector = (mphi & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
1838  fb = (mphi & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
1839  sectorFB = sector if fb == 0 else sector + 8
1840  tphi = phiTimes[phiKey]
1841  tphiTrunc = int(tphi) & 0x3ff
1842  for zKey in zTimes[layer]:
1843  mz = zKey & self.BKLM_MODULEID_MASK
1844  if mphi == mz:
1845  tz = zTimes[layer][zKey]
1846  tzTrunc = int(tz) & 0x3ff
1847  dt = (tphiTrunc - tzTrunc) & 0x3ff
1848  if dt >= 0x200:
1849  dt -= 0x400
1850  t = (tphi + tz) * 0.5
1851  tTrunc = int(t) & 0x3ff
1852  if layer < 2:
1853  self.hist_dtScint1d.Fill(dt)
1854  else:
1855  self.hist_dtRPC1d.Fill(dt)
1856  if abs(dt) < 4000:
1857  if layer < 2:
1858  self.hist_ctScintCal1d.Fill(tTrunc)
1859  else:
1860  self.hist_tRPCCal1d.Fill(tTrunc)
1861 
1862  # After processing all of the BKLMHit1ds in the event, draw the event display (perhaps)
1863 
1864  if (self.view == 1) and (self.eventDisplays < self.maxDisplays):
1865  drawnSectors = 0
1866  jCanvas = 1
1867  for sector in range(0, 8):
1868  if len(zyList[sector]) > self.minRPCHits:
1869  drawnSectors += 1
1870  self.eventCanvas.cd(jCanvas)
1871  title = 'e{0:02d}r{1}: event {2} z-readout hits in S{3}'.format(int(self.exp), int(self.run), event, sector)
1872  self.hist_ZY1D[jCanvas - 1].SetTitle(title)
1873  self.hist_ZY1D[jCanvas - 1].Draw()
1874  for g in self.bklmZY:
1875  g.Draw("L")
1876  for g in zyList[sector]:
1877  g.Draw("L")
1878  jCanvas += 1
1879  if jCanvas > 2:
1880  jCanvas = 1
1881  self.lastTitle = "Title:E{0} (#{1})".format(event, self.eventCounter)
1882  self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
1883  enoughXYHits = False
1884  for sector in range(0, 8):
1885  if len(xyList[sector]) > self.minRPCHits:
1886  enoughXYHits = True
1887  break
1888  if enoughXYHits:
1889  drawnSectors += 1
1890  self.eventCanvas.cd(jCanvas)
1891  jCanvas += 1
1892  title = 'e{0:02d}r{1}: event {2} phi-readout hits'.format(int(self.exp), int(self.run), event)
1893  self.hist_XY.SetTitle(title)
1894  self.hist_XY.Draw()
1895  for g in self.bklmXY:
1896  g.Draw("L")
1897  for sector in range(0, 8):
1898  for g in xyList[sector]:
1899  g.Draw("L")
1900  if jCanvas > 2:
1901  jCanvas = 1
1902  self.lastTitle = "Title:E{0} (#{1})".format(event, self.eventCounter)
1903  self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
1904  if jCanvas == 2:
1905  self.eventCanvas.cd(jCanvas)
1906  ROOT.gPad.Clear()
1907  self.lastTitle = "Title:E{0} (#{1})".format(event, self.eventCounter)
1908  self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
1909  if drawnSectors > 0:
1910  self.eventDisplays += 1
1911 
1912  # Process the BKLMHit2ds
1913 
1914  xyList = []
1915  zyList = []
1916  rpcHits = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
1917  for hit2d in hit2ds:
1918  key = hit2d.getModuleID()
1919  layer = (key & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
1920  sector = (key & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
1921  fb = (key & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
1922  phiStripMin = hit2d.getPhiStripMin() - 1
1923  phiStripMax = hit2d.getPhiStripMax() - 1
1924  zStripMin = hit2d.getZStripMin() - 1
1925  zStripMax = hit2d.getZStripMax() - 1
1926  sectorFB = sector if fb == 0 else sector + 8
1927  if layer >= 2:
1928  rpcHits[sectorFB] += 1
1929  if self.legacyTimes:
1930  dc = self.sectorFBToDC[sectorFB]
1931  copper = dc & 0x03
1932  finesse = dc >> 2
1933  n = rawklms[copper].GetDetectorNwords(0, finesse) >> 1
1934  trigCtime = (rawklms[copper].GetTTCtime(0) & 0x07ffffff) << 3
1935  ctDiffMax = 99999
1936  tCal = -1
1937  jZ = -1
1938  jPhi = -1
1939  ctZ = 0
1940  ctPhi = 0
1941  for j in range(0, n):
1942  if layer != rawLayer[dc][j]:
1943  continue
1944  if sector != rawSector[dc][j]:
1945  continue
1946  if fb != rawFb[dc][j]:
1947  continue
1948  strip = rawStrip[dc][j]
1949  plane = rawPlane[dc][j]
1950  if plane == 0: # it's a z strip
1951  if strip < zStripMin:
1952  continue
1953  if strip > zStripMax:
1954  continue
1955  ctZ = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
1956  jZ = j
1957  else: # it's a phi strip
1958  if strip < phiStripMin:
1959  continue
1960  if strip > phiStripMax:
1961  continue
1962  ctPhi = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
1963  jPhi = j
1964  if (jZ >= 0) and (jPhi >= 0):
1965  if layer < 2: # it's a scint layer
1966  if abs(ctZ - ctPhi) > 40:
1967  continue
1968  ct = int((ctZ + ctPhi) * 0.5 - trigCtime - self.ct0Scint[plane][sectorFB][layer]) & 0x3ff
1969  if abs(ct - self.ct0Cal) < ctDiffMax:
1970  ctDiffMax = int(abs(ct - self.ct0Cal))
1971  tCal = ct
1972  if ctDiffMax == 0:
1973  break
1974  else: # it's an RPC layer
1975  tCal = hit2d.getTime() - trigCtime - self.t0RPC[plane][sectorFB][layer]
1976  break
1977  else:
1978  if layer < 2:
1979  tCal = hit2d.getTime() - self.ct0Scint[plane][sectorFB][layer]
1980  else:
1981  tCal = hit2d.getTime() - self.t0RPC[plane][sectorFB][layer]
1982  tCalTrunc = int(tCal) & 0x3ff
1983  x = hit2d.getGlobalPositionX()
1984  y = hit2d.getGlobalPositionY()
1985  z = hit2d.getGlobalPositionZ()
1986  r = math.sqrt(x * x + y * y)
1987  isPromptHit = False
1988  promptColor = 3
1989  bkgdColor = 2
1990  if layer < 2:
1991  promptColor = 7
1992  bkgdColor = 4
1993  self.hist_ctScintCal2d.Fill(tCalTrunc)
1994  self.hist_ctScintCal2dBySector.Fill(sectorFB, tCalTrunc)
1995  if abs(tCalTrunc - self.ct0Cal2d) < 50:
1996  isPromptHit = True
1997  if fb == 0: # backward
1998  self.hist_occupancyBackwardXYPrompt.Fill(x, y)
1999  else: # forward
2000  self.hist_occupancyForwardXYPrompt.Fill(x, y)
2001  else:
2002  if fb == 0: # backward
2003  self.hist_occupancyBackwardXYBkgd.Fill(x, y)
2004  else: # forward
2005  self.hist_occupancyForwardXYBkgd.Fill(x, y)
2006  else:
2007  self.hist_tRPCCal2d.Fill(tCalTrunc)
2008  self.hist_tRPCCal2dBySector.Fill(sectorFB, tCalTrunc)
2009  if abs(tCalTrunc - self.t0Cal2d) < 50:
2010  isPromptHit = True
2011  self.hist_occupancyRZPrompt.Fill(z, layer)
2012  self.hist_occupancyZPrompt.Fill(z)
2013  self.hist_occupancyRPrompt.Fill(layer)
2014  if fb == 0: # backward
2015  self.hist_occupancyBackwardXYPrompt.Fill(x, y)
2016  self.hist_tVsZBwd.Fill(-(z - 47.0), tCalTrunc)
2017  else: # forward
2018  self.hist_occupancyForwardXYPrompt.Fill(x, y)
2019  self.hist_tVsZFwd.Fill(+(z - 47.0), tCalTrunc)
2020  elif abs(tCalTrunc - self.t0Cal2d) >= 50:
2021  self.hist_occupancyRZBkgd.Fill(z, layer)
2022  self.hist_occupancyZBkgd.Fill(z)
2023  self.hist_occupancyRBkgd.Fill(layer)
2024  if fb == 0: # backward
2025  self.hist_occupancyBackwardXYBkgd.Fill(x, y)
2026  else: # forward
2027  self.hist_occupancyForwardXYBkgd.Fill(x, y)
2028 
2029  # Add the hit to the event-display TGraph list (perhaps)
2030  if (self.view == 2) and (self.eventDisplays < self.maxDisplays):
2031  gXY = ROOT.TGraph()
2032  gXY.SetPoint(0, x - 1.0, y - 1.0)
2033  gXY.SetPoint(1, x - 1.0, y + 1.0)
2034  gXY.SetPoint(2, x + 1.0, y + 1.0)
2035  gXY.SetPoint(3, x + 1.0, y - 1.0)
2036  gXY.SetPoint(4, x - 1.0, y - 1.0)
2037  gXY.SetLineWidth(1)
2038  gZY = ROOT.TGraph()
2039  gZY.SetPoint(0, z - 1.0, y - 1.0)
2040  gZY.SetPoint(1, z - 1.0, y + 1.0)
2041  gZY.SetPoint(2, z + 1.0, y + 1.0)
2042  gZY.SetPoint(3, z + 1.0, y - 1.0)
2043  gZY.SetPoint(4, z - 1.0, y - 1.0)
2044  gZY.SetLineWidth(1)
2045  if isPromptHit:
2046  gXY.SetLineColor(promptColor)
2047  gZY.SetLineColor(promptColor)
2048  else:
2049  gXY.SetLineColor(bkgdColor)
2050  gZY.SetLineColor(bkgdColor)
2051  xyList.append(gXY)
2052  zyList.append(gZY)
2053 
2054  # After processing all of the hits in the event, draw the event display (perhaps)
2055 
2056  if (self.view == 2) and (self.eventDisplays < self.maxDisplays):
2057  hasEnoughRPCHits = False
2058  for count in rpcHits:
2059  if count > self.minRPCHits:
2060  hasEnoughRPCHits = True
2061  break
2062  if hasEnoughRPCHits:
2063  self.eventDisplays += 1
2064  title = 'e{0:02d}r{1}: event {2} 2D hits'.format(int(self.exp), int(self.run), event)
2065  self.hist_XY.SetTitle(title)
2066  self.hist_ZY.SetTitle(title)
2067  self.eventCanvas.cd(1)
2068  self.hist_XY.Draw()
2069  for g in self.bklmXY:
2070  g.Draw("L")
2071  for g in xyList:
2072  g.Draw("L")
2073  self.eventCanvas.cd(2)
2074  self.hist_ZY.Draw()
2075  for g in self.bklmZY:
2076  g.Draw("L")
2077  for g in zyList:
2078  g.Draw("L")
2079  self.lastTitle = "Title:E{0} (#{1})".format(event, self.eventCounter)
2080  self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
EventInspector.EventInspector.hist_unmappedChannelOccupancy
hist_unmappedChannelOccupancy
scatterplots of unmapped channel occupancy (1 hit per readout channel), indexed by sector#
Definition: EventInspector.py:378
EventInspector.EventInspector.BKLM_SECTOR_MASK
tuple BKLM_SECTOR_MASK
bit mask for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
Definition: EventInspector.py:44
EventInspector.EventInspector.hist_rawKLMchannelMultiplicity
hist_rawKLMchannelMultiplicity
scatterplots of multiplicity of entries in one readout channel vs lane/axis, indexed by sector#
Definition: EventInspector.py:278
EventInspector.EventInspector.hist_mappedScintTDCBySector
hist_mappedScintTDCBySector
scatterplot of scint mapped-channel TDC value (NOT relative to event's trigger Ctime) vs sector
Definition: EventInspector.py:536
EventInspector.EventInspector.BKLM_LAYER_MASK
tuple BKLM_LAYER_MASK
bit mask for layer-1 [0..15]; 0 is innermost and 14 is outermost
Definition: EventInspector.py:42
EventInspector.EventInspector.hist_XY
hist_XY
blank scatterplot to define the bounds of the BKLM end view
Definition: EventInspector.py:134
EventInspector.EventInspector.hist_rawKLMsizeByDCMultihit
hist_rawKLMsizeByDCMultihit
histograms of number of hits, including multiple entries on one readout channel, indexed by sector#
Definition: EventInspector.py:264
EventInspector.EventInspector.hist_tVsZBwd
hist_tVsZBwd
profile histogram of BKLMHit2d RPC time vs z (forward)
Definition: EventInspector.py:739
EventInspector.EventInspector.bklmXY
bklmXY
list of line-segment (x,y) points for the BKLM end view
Definition: EventInspector.py:757
EventInspector.EventInspector.hist_occupancyBackwardXYBkgd
hist_occupancyBackwardXYBkgd
scatterplot of end view of backward BKLM for out-of-time BKLMHit2ds
Definition: EventInspector.py:682
EventInspector.EventInspector.BKLM_MODULEID_MASK
tuple BKLM_MODULEID_MASK
bit mask for unique module identifier (end, sector, layer)
Definition: EventInspector.py:50
EventInspector.EventInspector.hist_mappedScintTime
hist_mappedScintTime
histogram of scint mapped-channel TDC value relative to event's trigger Ctime
Definition: EventInspector.py:532
EventInspector.EventInspector.hist_occupancyZBkgd
hist_occupancyZBkgd
histogram of z coordinate for out-of-time BKLMHit2ds
Definition: EventInspector.py:702
EventInspector.EventInspector.hist_occupancyRPrompt
hist_occupancyRPrompt
histogram of layer# for in-time BKLMHit2ds
Definition: EventInspector.py:694
EventInspector.EventInspector.hist_tphiRPCCal1d
hist_tphiRPCCal1d
histogram of RPC-phi BKLMHit1d time relative to event's trigger time, corrected for inter-sector vari...
Definition: EventInspector.py:625
EventInspector.EventInspector.hist_mappedSectorOccupancyMultihit
hist_mappedSectorOccupancyMultihit
histogram of number of mapped hits by sector, including multiple entries on one readout channel
Definition: EventInspector.py:294
EventInspector.EventInspector.hist_ctimeRangeVsTimeRPC
hist_ctimeRangeVsTimeRPC
scatterplot of RPC Ctime range vs time
Definition: EventInspector.py:588
EventInspector.EventInspector.ct0Cal1d
ct0Cal1d
scint-ctime calibration adjustment (ns) for BKLMHit1ds
Definition: EventInspector.py:164
EventInspector.EventInspector.hist_nHit1dRPCPrompt
hist_nHit1dRPCPrompt
histogram of the number of in-time RPC BKLMHit1ds
Definition: EventInspector.py:599
EventInspector.EventInspector.hist_rawKLMtdcExtraScint
hist_rawKLMtdcExtraScint
scatterplot of the RawKLM scint hit's extra bits vs sector in the third (time) word
Definition: EventInspector.py:248
EventInspector.EventInspector.endRun
def endRun(self)
Definition: EventInspector.py:1409
EventInspector.EventInspector.hist_ctScintCal2d
hist_ctScintCal2d
histogram of scint calibrated time in BKLMHit2ds
Definition: EventInspector.py:721
EventInspector.EventInspector.hist_rawKLMchannelMultiplicityFine
hist_rawKLMchannelMultiplicityFine
scatterplots of multiplicity of entries in one readout channel vs lane/axis/channel,...
Definition: EventInspector.py:280
EventInspector.EventInspector.hist_unmappedScintTimeBySector
hist_unmappedScintTimeBySector
scatterplot of scint unmapped-channel TDC value relative to event's trigger Ctime vs sector
Definition: EventInspector.py:553
EventInspector.EventInspector.hist_unmappedScintCtime
hist_unmappedScintCtime
histogram of scint unmapped-channel CTIME value relative to event's trigger Ctime
Definition: EventInspector.py:516
EventInspector.EventInspector.hist_mappedScintSectorOccupancy
hist_mappedScintSectorOccupancy
scatterplot of number of mapped scint hits by lane/axis vs sector, at most one entry per readout chan...
Definition: EventInspector.py:344
EventInspector.EventInspector.hist_occupancyForwardXYPrompt
hist_occupancyForwardXYPrompt
scatterplot of end view of forward BKLM for in-time BKLMHit2ds
Definition: EventInspector.py:670
EventInspector.EventInspector.hist_ctimeRangeRPC
hist_ctimeRangeRPC
histogram of RPC Ctime range
Definition: EventInspector.py:571
EventInspector.EventInspector.hist_mappedRPCCtimeRange
hist_mappedRPCCtimeRange
histogram of RPC mapped-channel REVO9 range in event
Definition: EventInspector.py:439
EventInspector.EventInspector.hist_mappedScintCtime
hist_mappedScintCtime
histogram of scint mapped-channel CTIME value relative to event's trigger Ctime
Definition: EventInspector.py:467
EventInspector.EventInspector.maxDisplays
maxDisplays
internal copy of the maximum number of event displays to write
Definition: EventInspector.py:82
EventInspector.EventInspector.hist_dtRPC1d
hist_dtRPC1d
histogram of RPC-phi and -z BKLMHit1d time difference
Definition: EventInspector.py:640
EventInspector.EventInspector.beginRun
def beginRun(self)
Definition: EventInspector.py:1404
EventInspector.EventInspector.hist_nHit1d
hist_nHit1d
histogram of the number of BKLMHit1ds
Definition: EventInspector.py:597
EventInspector.EventInspector.BKLM_ID
int BKLM_ID
COPPER base identifier for BKLM readout.
Definition: EventInspector.py:22
EventInspector.EventInspector.t0Cal
t0Cal
Time-calibration constants obtained from experiment 7 run 1505 RPC-time calibration adjustment (ns) f...
Definition: EventInspector.py:156
EventInspector.EventInspector.hist_occupancyForwardXYBkgd
hist_occupancyForwardXYBkgd
scatterplot of end view of forward BKLM for out-of-time BKLMHit2ds
Definition: EventInspector.py:678
EventInspector.EventInspector.hist_mappedRPCCtimeRangeBySector
hist_mappedRPCCtimeRangeBySector
scatterplot of RPC mapped-channel REVO9 range in event vs sector
Definition: EventInspector.py:444
EventInspector.EventInspector.t0Cal1d
t0Cal1d
RPC-time calibration adjustment (ns) for BKLMHit1ds.
Definition: EventInspector.py:158
EventInspector.EventInspector.BKLM_SECTION_MASK
tuple BKLM_SECTION_MASK
bit mask for section [0..1]; forward is 0
Definition: EventInspector.py:46
EventInspector.EventInspector.hist_tRPCCal1d
hist_tRPCCal1d
histogram of RPC-phi and -z BKLMHit1d avg time relative to event's trigger time, corrected for inter-...
Definition: EventInspector.py:635
EventInspector.EventInspector.BKLM_SECTION_BIT
int BKLM_SECTION_BIT
bit position for section [0..1]; forward is 0
Definition: EventInspector.py:34
EventInspector.EventInspector.hist_n1dPhiZ
hist_n1dPhiZ
scatterplot of #Z BKLMHit1ds vs #Phi BKLMHit1ds
Definition: EventInspector.py:609
EventInspector.EventInspector.hist_rawKLMnumNodes
hist_rawKLMnumNodes
histogram of the RawKLM's NumNodes (should be 1)
Definition: EventInspector.py:222
Belle2::PyStoreObj
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:69
EventInspector.EventInspector.hist_unmappedRPCTimeBySector
hist_unmappedRPCTimeBySector
scatterplot of RPC unmapped-channel TDC value relative to event's trigger time, by sector
Definition: EventInspector.py:455
EventInspector.EventInspector.hist_tRPCCal2d
hist_tRPCCal2d
histogram of RPC calibrated time in BKLMHit2ds
Definition: EventInspector.py:710
EventInspector.EventInspector.hist_unmappedScintCtimeBySector
hist_unmappedScintCtimeBySector
scatterplot of scint unmapped-channel CTIME value relative to event's trigger Ctime,...
Definition: EventInspector.py:521
EventInspector.EventInspector.hist_mappedScintCtimeCal
hist_mappedScintCtimeCal
histogram of scint mapped-channel CTIME value relative to event's trigger Ctime, corrected for inter-...
Definition: EventInspector.py:478
EventInspector.EventInspector.EKLM_ID
int EKLM_ID
COPPER base identifier for EKLM readout.
Definition: EventInspector.py:24
EventInspector.EventInspector.BKLM_LAYER_BIT
int BKLM_LAYER_BIT
bit position for layer-1 [0..14]; 0 is innermost
Definition: EventInspector.py:30
EventInspector.EventInspector.hist_RPCTimeLowBitsBySector
hist_RPCTimeLowBitsBySector
scatterplot of RPC TDC low-order bits vs sector (should be 0 since granularity is 4 ns)
Definition: EventInspector.py:401
EventInspector.EventInspector.hist_ZY
hist_ZY
blank scatterplot to define the bounds of the BKLM side view for 2D hits
Definition: EventInspector.py:143
EventInspector.EventInspector.hist_rawKLMsizeMultihit
hist_rawKLMsizeMultihit
histogram of number of hits, including multiple entries on one readout channel
Definition: EventInspector.py:260
EventInspector.EventInspector.hist_trigCtimeVsTrigRevo9time
hist_trigCtimeVsTrigRevo9time
histogram of RawKLM[] header's trigger CTIME relative to its final-data-word trigger REVO9 time
Definition: EventInspector.py:562
EventInspector.EventInspector.hist_mappedRPCTime
hist_mappedRPCTime
histogram of RPC mapped-channel TDC value relative to event's trigger time
Definition: EventInspector.py:407
EventInspector.EventInspector.hist_mappedScintLaneAxisOccupancy
hist_mappedScintLaneAxisOccupancy
scatterplot of number of mapped scint hits by lane/axis vs sector, at most one entry per readout chan...
Definition: EventInspector.py:350
EventInspector.EventInspector.hist_occupancyRBkgd
hist_occupancyRBkgd
histogram of layer# for out-of-time BKLMHit2ds
Definition: EventInspector.py:706
EventInspector.EventInspector.hist_rawKLMnodeID
hist_rawKLMnodeID
scatterplot of the RawKLM's COPPER index vs NodeID relative to the base BKLM/EKLM values
Definition: EventInspector.py:224
EventInspector.EventInspector.hist_mappedRPCSectorOccupancy
hist_mappedRPCSectorOccupancy
histogram of number of mapped RPC hits by sector, at most one entry per readout channel
Definition: EventInspector.py:318
EventInspector.EventInspector.pdfName
pdfName
internal copy of the pathname of the output histogram PDF file
Definition: EventInspector.py:76
EventInspector.EventInspector.hist_multiplicityZBySector
hist_multiplicityZBySector
scatterplot of #Z BKLMHit1ds vs sector
Definition: EventInspector.py:619
EventInspector.EventInspector.BKLM_SECTOR_BIT
int BKLM_SECTOR_BIT
bit position for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
Definition: EventInspector.py:32
EventInspector.EventInspector.hist_ScintTimeLowBitsBySector
hist_ScintTimeLowBitsBySector
scatterplot of scint TDC low-order bits vs sector
Definition: EventInspector.py:461
EventInspector.EventInspector.makeText
def makeText(self, x, y, s)
Definition: EventInspector.py:113
EventInspector.EventInspector.eventCounter
eventCounter
event counter (needed for PDF table of contents' ordinal event#)
Definition: EventInspector.py:92
EventInspector.EventInspector.hist_mappedChannelOccupancyPrompt
hist_mappedChannelOccupancyPrompt
scatterplots of in-time mapped channel occupancy (1 hit per readout channel), indexed by sector#
Definition: EventInspector.py:370
EventInspector.EventInspector.hist_mappedChannelOccupancyBkgd
hist_mappedChannelOccupancyBkgd
scatterplots of out-of-time mapped channel occupancy (1 hit per readout channel), indexed by sector#
Definition: EventInspector.py:374
EventInspector.EventInspector.hist_ctScintCal1d
hist_ctScintCal1d
histogram of scint-phi and -z BKLMHit1d avg time relative to event's trigger Ctime,...
Definition: EventInspector.py:655
EventInspector.EventInspector.event
def event(self)
Definition: EventInspector.py:1414
EventInspector.EventInspector.initialize
def initialize(self)
Definition: EventInspector.py:129
EventInspector.EventInspector.ct0Cal2d
ct0Cal2d
scint-ctime calibration adjustment (ns) for BKLMHit2ds
Definition: EventInspector.py:166
EventInspector.EventInspector.run
run
internal copy of run number
Definition: EventInspector.py:72
EventInspector.EventInspector.hist_nHit1dPrompt
hist_nHit1dPrompt
histogram of the number of in-time scint BKLMHit1ds
Definition: EventInspector.py:605
EventInspector.EventInspector.hist_occupancyBackwardXYPrompt
hist_occupancyBackwardXYPrompt
scatterplot of end view of backward BKLM for in-time BKLMHit2ds
Definition: EventInspector.py:674
EventInspector.EventInspector.hist_tdcRangeRPC
hist_tdcRangeRPC
histogram of RPC TDC range
Definition: EventInspector.py:566
EventInspector.EventInspector.hist_mappedRPCTimeCal
hist_mappedRPCTimeCal
histogram of RPC mapped-channel TDC value relative to event's trigger time, corrected for inter-secto...
Definition: EventInspector.py:410
EventInspector.EventInspector.hist_mappedScintZCtimePerLayer
hist_mappedScintZCtimePerLayer
histograms of scint mapped-channel z-strip CTIME value relative to event's trigger Ctime,...
Definition: EventInspector.py:491
EventInspector.EventInspector.BKLM_STRIP_MASK
int BKLM_STRIP_MASK
bit mask for strip-1 [0..47]
Definition: EventInspector.py:38
EventInspector.EventInspector.hist_nHit2d
hist_nHit2d
histogram of the number of BKLMHit2ds
Definition: EventInspector.py:668
bklmDB.fillDB
def fillDB()
Definition: bklmDB.py:9
EventInspector.EventInspector.hist_mappedRPCLaneAxisOccupancy
hist_mappedRPCLaneAxisOccupancy
scatterplot of number of mapped RPC hits by lane/axis vs sector, at most one entry per readout channe...
Definition: EventInspector.py:324
EventInspector.EventInspector.hist_ctzScintCal1d
hist_ctzScintCal1d
histogram of scint-z BKLMHit1d time relative to event's trigger Ctime, corrected for inter-sector var...
Definition: EventInspector.py:650
EventInspector.EventInspector.hist_mappedSectorOccupancy
hist_mappedSectorOccupancy
histogram of number of mapped hits by sector, at most one entry per readout channel
Definition: EventInspector.py:306
EventInspector.EventInspector.dcToSectorFB
dcToSectorFB
map for data concentrator -> sectorFB
Definition: EventInspector.py:153
EventInspector.EventInspector.bklmZY
bklmZY
list of line-segment (z,y) points for the BKLM side view
Definition: EventInspector.py:806
EventInspector.EventInspector.hist_mappedRPCPhiTimePerLayer
hist_mappedRPCPhiTimePerLayer
histograms of RPC mapped-channel phi-strip TDC value relative to event's trigger time,...
Definition: EventInspector.py:413
EventInspector.EventInspector.view
view
view event displays using one-dimensional (1) or two-dimensional (2) hits
Definition: EventInspector.py:90
EventInspector.EventInspector.t0RPC
t0RPC
per-layer variations in RPC z- and phi-time calibration adjustment (ns) for rawKLMs
Definition: EventInspector.py:168
EventInspector.EventInspector.makeGraph
def makeGraph(self, x, y)
Definition: EventInspector.py:99
EventInspector.EventInspector.hist_rawKLMadcExtraScint
hist_rawKLMadcExtraScint
scatterplot of the RawKLM scint hit's extra bits vs sector in the fourth (adc) word
Definition: EventInspector.py:254
EventInspector.EventInspector.exp
exp
internal copy of experiment number
Definition: EventInspector.py:70
EventInspector.EventInspector.BKLM_STRIP_BIT
int BKLM_STRIP_BIT
bit position for strip-1 [0..47]
Definition: EventInspector.py:26
EventInspector.EventInspector.hist_tRPCCal2dBySector
hist_tRPCCal2dBySector
scatterplot of RPC calibrated time in BKLMHit2ds vs sector
Definition: EventInspector.py:715
EventInspector.EventInspector.ct0Cal
ct0Cal
scint-ctime calibration adjustment (ns) for rawKLMs
Definition: EventInspector.py:162
EventInspector.EventInspector.hist_unmappedScintTime
hist_unmappedScintTime
histogram of scint unmapped-channel TDC value relative to event's trigger Ctime
Definition: EventInspector.py:548
EventInspector.EventInspector.verbosity
verbosity
internal copy of the histogram verbosity in the histogram PDF file
Definition: EventInspector.py:80
EventInspector.EventInspector.hist_multiplicityPhiBySector
hist_multiplicityPhiBySector
scatterplot of #Phi BKLMHit1ds vs sector
Definition: EventInspector.py:613
EventInspector.EventInspector.BKLM_MAXSTRIP_BIT
int BKLM_MAXSTRIP_BIT
bit position for maxStrip-1 [0..47]
Definition: EventInspector.py:36
EventInspector.EventInspector.hist_tzRPCCal1d
hist_tzRPCCal1d
histogram of RPC-z BKLMHit1d time relative to event's trigger time, corrected for inter-sector variat...
Definition: EventInspector.py:630
EventInspector.EventInspector.hist_tdcRangeVsCtimeRangeRPC
hist_tdcRangeVsCtimeRangeRPC
scatterplot of RPC TDC range vs Ctime range
Definition: EventInspector.py:576
EventInspector.EventInspector.electIdToModuleId
electIdToModuleId
readout <-> detector map (from the information retrieved from the conditions database)
Definition: EventInspector.py:149
EventInspector.EventInspector.singleEntry
singleEntry
select events with any (0) or exactly one (1) or more than one (2) entries/channel
Definition: EventInspector.py:88
EventInspector.EventInspector.hist_unmappedSectorOccupancyMultihit
hist_unmappedSectorOccupancyMultihit
histogram of number of unmapped hits by sector, including multiple entries on one readout channel
Definition: EventInspector.py:300
EventInspector.EventInspector.hist_rawKLMsize
hist_rawKLMsize
histogram of number of hits, at most one entry per readout channel
Definition: EventInspector.py:262
EventInspector.EventInspector
Definition: EventInspector.py:17
EventInspector.EventInspector.hist_mappedRPCTimeCalBySector
hist_mappedRPCTimeCalBySector
scatterplot of RPC mapped-channel TDC relative to trigger time, corrected for inter-sector variation,...
Definition: EventInspector.py:433
EventInspector.EventInspector.legacyTimes
legacyTimes
calculate prompt time for legacy BKLMHit1ds and BKLMHit2ds (True) or use stored time (False)
Definition: EventInspector.py:86
EventInspector.EventInspector.hist_mappedScintPhiCtimePerLayer
hist_mappedScintPhiCtimePerLayer
histograms of scint mapped-channel phi-strip CTIME value relative to event's trigger Ctime,...
Definition: EventInspector.py:489
EventInspector.EventInspector.hist_occupancyZPrompt
hist_occupancyZPrompt
histogram of z coordinate for in-time BKLMHit2ds
Definition: EventInspector.py:690
EventInspector.EventInspector.hist_ctScintCal2dBySector
hist_ctScintCal2dBySector
scatterplot of scint calibrated time in BKLMHit2ds vs sector
Definition: EventInspector.py:726
EventInspector.EventInspector.__init__
def __init__(self, exp, run, histName, pdfName, eventPdfName, verbosity, maxDisplays, minRPCHits, legacyTimes, singleEntry, view)
Definition: EventInspector.py:52
EventInspector.EventInspector.hist_mappedScintTDC
hist_mappedScintTDC
histogram of scint mapped-channel TDC value (NOT relative to event's trigger Ctime)
Definition: EventInspector.py:527
EventInspector.EventInspector.hist_unmappedScintSectorOccupancy
hist_unmappedScintSectorOccupancy
histogram of number of unmapped scint hits by sector, at most one entry per readout channel
Definition: EventInspector.py:357
EventInspector.EventInspector.eventCanvas
eventCanvas
TCanvas on which event displays will be drawn.
Definition: EventInspector.py:748
EventInspector.EventInspector.hist_rawKLMadcExtraRPC
hist_rawKLMadcExtraRPC
scatterplot of the RawKLM RPC hit's extra bits vs sector in the fourth (adc) word
Definition: EventInspector.py:242
EventInspector.EventInspector.hist_mappedScintCtimeRangeBySector
hist_mappedScintCtimeRangeBySector
scatterplot of scint mapped-channel CTIME range in event vs sector
Definition: EventInspector.py:510
EventInspector.EventInspector.hist_mappedScintCtimeCalBySector
hist_mappedScintCtimeCalBySector
scatterplot of scint mapped-channel CTIME relative to trigger Ctime, corrected for inter-sector varia...
Definition: EventInspector.py:483
EventInspector.EventInspector.hist_rawKLMlaneFlag
hist_rawKLMlaneFlag
scatterplot of the RawKLM hit's lane vs flag (1=RPC, 2=Scint)
Definition: EventInspector.py:230
EventInspector.EventInspector.hist_mappedRPCZTimePerLayer
hist_mappedRPCZTimePerLayer
histograms of RPC mapped-channel z-strip TDC value relative to event's trigger time,...
Definition: EventInspector.py:415
Belle2::PyStoreArray
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:58
EventInspector.EventInspector.eventPdfName
eventPdfName
internal copy of the pathname of the output event-display PDF file
Definition: EventInspector.py:78
EventInspector.EventInspector.hist_nHit1dScint
hist_nHit1dScint
histogram of the number of scint BKLMHit1ds
Definition: EventInspector.py:603
EventInspector.EventInspector.t0Cal2d
t0Cal2d
RPC-time calibration adjustment (ns) for BKLMHit2ds.
Definition: EventInspector.py:160
EventInspector.EventInspector.hist_occupancyRZBkgd
hist_occupancyRZBkgd
scatterplot of side view of forward BKLM for in-time BKLMHit2ds
Definition: EventInspector.py:698
EventInspector.EventInspector.ct0Scint
ct0Scint
per-layer variations in scint-ctime calibration adjustment (ns) for rawKLMs
Definition: EventInspector.py:202
EventInspector.EventInspector.hist_ctphiScintCal1d
hist_ctphiScintCal1d
histogram of scint-phi BKLMHit1d time relative to event's trigger Ctime, corrected for inter-sector v...
Definition: EventInspector.py:645
EventInspector.EventInspector.hist_unmappedRPCSectorOccupancy
hist_unmappedRPCSectorOccupancy
histogram of number of unmapped RPC hits by sector, at most one entry per readout channel
Definition: EventInspector.py:331
EventInspector.EventInspector.hist_mappedScintCtimeRange
hist_mappedScintCtimeRange
histogram of scint mapped-channel CTIME range in event
Definition: EventInspector.py:505
EventInspector.EventInspector.hist_mappedScintTimeBySector
hist_mappedScintTimeBySector
scatterplot of scint mapped-channel TDC value relative to event's trigger Ctime vs sector
Definition: EventInspector.py:542
EventInspector.EventInspector.hist_ZY1D
hist_ZY1D
blank scatterplot to define the bounds of the BKLM side view for 1D hits
Definition: EventInspector.py:137
EventInspector.EventInspector.hist_unmappedRPCTime
hist_unmappedRPCTime
histogram of RPC unmapped-channel TDC value relative to event's trigger time
Definition: EventInspector.py:450
EventInspector.EventInspector.hist_dtScint1d
hist_dtScint1d
histogram of scint-phi and -z BKLMHit1d time difference
Definition: EventInspector.py:660
EventInspector.EventInspector.hist_nDigit
hist_nDigit
histogram of the number of BKLMDigits in the event
Definition: EventInspector.py:216
EventInspector.EventInspector.hist_nRawKLM
hist_nRawKLM
histogram of the number of RawKLMs in the event (should be 1)
Definition: EventInspector.py:218
EventInspector.EventInspector.hist_nHit1dBkgd
hist_nHit1dBkgd
histogram of the number of out-of-time scint BKLMHit1ds
Definition: EventInspector.py:607
EventInspector.EventInspector.sectorFBToDC
sectorFBToDC
map for sectorFB -> data concentrator
Definition: EventInspector.py:151
EventInspector.EventInspector.terminate
def terminate(self)
Definition: EventInspector.py:890
EventInspector.EventInspector.hist_unmappedSectorOccupancy
hist_unmappedSectorOccupancy
histogram of number of unmapped hits by sector, at most one entry per readout channel
Definition: EventInspector.py:312
EventInspector.EventInspector.hist_mappedScintCtimeBySector
hist_mappedScintCtimeBySector
scatterplot of scint mapped-channel CTIME value relative to event's trigger Ctime vs sector
Definition: EventInspector.py:472
EventInspector.EventInspector.hist_rawKLMnumEvents
hist_rawKLMnumEvents
histogram of the RawKLM's NumEvents (should be 1)
Definition: EventInspector.py:220
EventInspector.EventInspector.eventDisplays
eventDisplays
event-display counter
Definition: EventInspector.py:94
EventInspector.EventInspector.histogramFile
histogramFile
Output ROOT TFile that will contain the histograms/scatterplots.
Definition: EventInspector.py:208
EventInspector.EventInspector.minRPCHits
minRPCHits
internal copy of the minimum number of RPC BKLMHit2ds in any sector for event display
Definition: EventInspector.py:84
EventInspector.EventInspector.hist_tdcRangeVsTimeRPC
hist_tdcRangeVsTimeRPC
scatterplot of RPC TDC range vs time
Definition: EventInspector.py:582
EventInspector.EventInspector.hist_rawKLMtdcExtraRPC
hist_rawKLMtdcExtraRPC
scatterplot of the RawKLM RPC hit's extra bits vs sector in the third (time) word
Definition: EventInspector.py:236
EventInspector.EventInspector.lastTitle
lastTitle
title of the last-drawn event display (needed for PDF table of contents' last event)
Definition: EventInspector.py:96
EventInspector.EventInspector.histName
histName
internal copy of the pathname of the output histogram ROOT file
Definition: EventInspector.py:74
EventInspector.EventInspector.hist_nHit1dRPCBkgd
hist_nHit1dRPCBkgd
histogram of the number of out-of-time RPC BKLMHit1ds
Definition: EventInspector.py:601
EventInspector.EventInspector.BKLM_PLANE_BIT
int BKLM_PLANE_BIT
bit position for plane-1 [0..1]; 0 is inner-plane
Definition: EventInspector.py:28
EventInspector.EventInspector.hist_unmappedRPCLaneAxisOccupancy
hist_unmappedRPCLaneAxisOccupancy
scatterplot of number of unmapped RPC hits by lane/axis vs sector, at most one entry per readout chan...
Definition: EventInspector.py:337
EventInspector.EventInspector.BKLM_PLANE_MASK
tuple BKLM_PLANE_MASK
bit mask for plane-1 [0..1]; 0 is inner-plane
Definition: EventInspector.py:40
EventInspector.EventInspector.hist_rawKLMsizeByDC
hist_rawKLMsizeByDC
histograms of number of hits, at most one entry per readout channel, indexed by sector#
Definition: EventInspector.py:266
EventInspector.EventInspector.hist_tVsZFwd
hist_tVsZFwd
profile histogram of BKLMHit2d RPC time vs z (forward)
Definition: EventInspector.py:733
EventInspector.EventInspector.BKLM_MAXSTRIP_MASK
tuple BKLM_MAXSTRIP_MASK
bit mask for maxStrip-1 [0..47]
Definition: EventInspector.py:48
EventInspector.EventInspector.hist_occupancyRZPrompt
hist_occupancyRZPrompt
scatterplot of side view of forward BKLM for in-time BKLMHit2ds
Definition: EventInspector.py:686
EventInspector.EventInspector.hist_unmappedScintLaneAxisOccupancy
hist_unmappedScintLaneAxisOccupancy
scatterplot of number of unmapped scint hits by lane/axis vs sector, at most one entry per readout ch...
Definition: EventInspector.py:363
EventInspector.EventInspector.hist_mappedRPCTimeBySector
hist_mappedRPCTimeBySector
scatterplot of RPC mapped-channel TDC value relative to event's trigger time vs sector
Definition: EventInspector.py:427