Belle II Software development
EventInspectorLookback Class Reference
Inheritance diagram for EventInspectorLookback:

Public Member Functions

def __init__ (self, exp, run, histName, pdfName, mode, window)
 
def initialize (self)
 
def terminate (self)
 
def beginRun (self)
 
def endRun (self)
 
def event (self)
 

Public Attributes

 exp
 internal copy of experiment number
 
 run
 internal copy of run number
 
 histName
 internal copy of the pathname of the output histogram ROOT file
 
 pdfName
 internal copy of the pathname of the output histogram PDF file
 
 windowMode
 window mode as a string for histogram labels/titles
 
 windowMinValue
 highest observed lookback-window value
 
 windowMaxValue
 highest observed lookback-window value
 
 windowStepValue
 lookback-window value step
 
 electIdToModuleId
 readout <-> detector map (from the information retrieved from the conditions database)
 
 sectorFBToDC
 map for sectorFB -> data concentrator
 
 dcToSectorFB
 map for data concentrator -> sectorFB
 
 histogramFile
 Output ROOT TFile that will contain the histograms/scatterplots.
 
 hist_mappedRPCTimeCal
 histogram of RPC TDC - trigger value
 
 dict_mappedRPCTimeCalByWindow
 dictionary of histograms of RPC TDC - trigger value, keyed by lookback-window value
 
 dict_nRawKLMs
 dictionary of the number of RawKLM hits for each lookback-window value
 
 hist_mappedRPCTimeCalByWindow
 reference to the RPC-time histogram for the currevent value of the lookback window parameter
 
 hist_occupancyForwardXY
 scatterplot of end view of forward BKLM for all KLMHit2ds
 
 hist_occupancyBackwardXY
 scatterplot of end view of backward BKLM for all KLMHit2ds
 
 hist_occupancyXYByWindow
 reference to the xy scatterplot for the currevent value of the lookback window parameter
 
 dict_occupancyXYByWindow
 dictionary of scatterplots of end view of forward BKLM, keyed by lookback-window value
 
 dict_nHit2ds
 dictionary of the number of KLMHit2ds for each lookback-window value
 
 dict_nEvents
 dictionary of the number of events for each lookback-window value, for normalization
 
 windowValue
 cached value of the lookback-window value, to avoid unnecessary reassignments-to-same-value in event()
 

Static Public Attributes

int BKLM_ID = 0x07000000
 COPPER base identifier for BKLM readout.
 
int EKLM_ID = 0x08000000
 COPPER base identifier for EKLM readout.
 
int BKLM_STRIP_BIT = 0
 bit position for strip-1 [0..47]
 
int BKLM_PLANE_BIT = 6
 bit position for plane-1 [0..1]; 0 is inner-plane
 
int BKLM_LAYER_BIT = 7
 bit position for layer-1 [0..14]; 0 is innermost
 
int BKLM_SECTOR_BIT = 11
 bit position for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
 
int BKLM_SECTION_BIT = 14
 bit position for section [0..1]; forward is 0
 
int BKLM_MAXSTRIP_BIT = 15
 bit position for maxStrip-1 [0..47]
 
int BKLM_STRIP_MASK = 0x3f
 bit mask for strip-1 [0..47]
 
tuple BKLM_PLANE_MASK = (1 << BKLM_PLANE_BIT)
 bit mask for plane-1 [0..1]; 0 is inner-plane
 
tuple BKLM_LAYER_MASK = (15 << BKLM_LAYER_BIT)
 bit mask for layer-1 [0..15]; 0 is innermost and 14 is outermost
 
tuple BKLM_SECTOR_MASK = (7 << BKLM_SECTOR_BIT)
 bit mask for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis
 
tuple BKLM_SECTION_MASK = (1 << BKLM_SECTION_BIT)
 bit mask for section [0..1]; forward is 0
 
tuple BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
 bit mask for maxStrip-1 [0..47]
 
tuple BKLM_MODULEID_MASK = (BKLM_SECTION_MASK | BKLM_SECTOR_MASK | BKLM_LAYER_MASK)
 bit mask for unique module identifier (end, sector, layer)
 

Detailed Description

Analyze RPC lookback-window parameter settings, fill histograms

Definition at line 23 of file EventInspectorLookback.py.

Constructor & Destructor Documentation

◆ __init__()

def __init__ (   self,
  exp,
  run,
  histName,
  pdfName,
  mode,
  window 
)
Constructor

Arguments:
    exp (str): formatted experiment number
    run (str): formatter run number
    histName (str): path name of the output histogram ROOT file
    pdfName (str): path name of the output histogram PDF file
    mode (int): specifies the lookback-window mode
                0: coarse window start values
                1: coarse window width values
                2: fine window start values
                3: fine window width values
    window (int, int, int): specifies the lookback-window min, max and step values

Definition at line 57 of file EventInspectorLookback.py.

57 def __init__(self, exp, run, histName, pdfName, mode, window):
58 """Constructor
59
60 Arguments:
61 exp (str): formatted experiment number
62 run (str): formatter run number
63 histName (str): path name of the output histogram ROOT file
64 pdfName (str): path name of the output histogram PDF file
65 mode (int): specifies the lookback-window mode
66 0: coarse window start values
67 1: coarse window width values
68 2: fine window start values
69 3: fine window width values
70 window (int, int, int): specifies the lookback-window min, max and step values
71 """
72 super().__init__()
73
74 self.exp = exp
75
76 self.run = run
77
78 self.histName = histName
79
80 self.pdfName = pdfName
81
82 windowModes = {0: "coarse start", 1: "coarse width", 2: "fine start", 3: "fine width"}
83
84 self.windowMode = windowModes[mode]
85
86 self.windowMinValue = window[0]
87
88 self.windowMaxValue = window[1]
89
90 self.windowStepValue = window[2]
91 print(f"Mode = {mode} start = {window[0]} end = {window[1]} step = {window[2]}")
92

Member Function Documentation

◆ beginRun()

def beginRun (   self)
Handle begin of run: print diagnostic message

Definition at line 252 of file EventInspectorLookback.py.

252 def beginRun(self):
253 """Handle begin of run: print diagnostic message"""
254 EventMetaData = Belle2.PyStoreObj('EventMetaData')
255 print('beginRun', EventMetaData.getRun())
256
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67

◆ endRun()

def endRun (   self)
Handle end of run: print diagnostic message

Definition at line 257 of file EventInspectorLookback.py.

257 def endRun(self):
258 """Handle end of run: print diagnostic message"""
259 EventMetaData = Belle2.PyStoreObj('EventMetaData')
260 print('endRun', EventMetaData.getRun())
261

◆ event()

def event (   self)
Process one event: fill histograms

Definition at line 262 of file EventInspectorLookback.py.

262 def event(self):
263 """Process one event: fill histograms"""
264
265 EventMetaData = Belle2.PyStoreObj('EventMetaData')
266 event = EventMetaData.getEvent()
267 rawklms = Belle2.PyStoreArray('RawKLMs')
268 hit2ds = Belle2.PyStoreArray('KLMHit2ds')
269
270 # Process the RawKLMs
271
272 for copper in range(0, len(rawklms)):
273 rawklm = rawklms[copper]
274 if rawklm.GetNumEntries() != 1:
275 print('##0 Event', event, 'copper', copper, ' getNumEntries=', rawklm.GetNumEntries())
276 continue
277 nodeID = rawklm.GetNodeID(0) - self.BKLM_ID
278 if nodeID >= self.EKLM_ID - self.BKLM_ID:
279 nodeID = nodeID - (self.EKLM_ID - self.BKLM_ID) + 4
280 if (nodeID < 0) or (nodeID > 4): # examine BKLM nodes only
281 continue
282 trigCtime = (rawklm.GetTTCtime(0) & 0x7ffffff) << 3 # (ns)
283 for finesse in range(0, 4):
284 dc = (finesse << 2) + (copper & 0x3)
285 sectorFB = self.dcToSectorFB[dc] # noqa (F841) kept for completeness
286 nWords = rawklm.GetDetectorNwords(0, finesse)
287 if nWords <= 0:
288 continue
289 bufSlot = rawklm.GetDetectorBuffer(0, finesse)
290 lastWord = bufSlot[nWords - 1]
291 windowValue = (lastWord >> 16) & 0xffff
292 if windowValue != self.windowValue:
293 if windowValue in self.dict_nEvents:
294 self.windowValue = windowValue
295 self.hist_mappedRPCTimeCalByWindow = self.dict_mappedRPCTimeCalByWindow[windowValue]
296 self.hist_occupancyXYByWindow = self.dict_occupancyXYByWindow[windowValue]
297 else:
298 return # skip bogus event, including event with seed value of 0xcafe
299 if lastWord & 0xffff != 0:
300 print("##1 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'lastWord=', hex(lastWord))
301 if (nWords % 2) == 0:
302 print("##2 Event", event, 'copper', copper, 'finesse', finesse, 'n=', nWords, 'should be odd -- skipping')
303 continue
304 n = nWords >> 1 # number of Data-Concentrator data packets
305 self.dict_nRawKLMs[self.windowValue] += n
306 # first (and only) pass over this DC's hits: histogram everything
307 for j in range(0, n):
308 word0 = bufSlot[j * 2]
309 word1 = bufSlot[j * 2 + 1]
310 ctime = word0 & 0xffff # noqa (F841) kept for completeness
311 channel = (word0 >> 16) & 0x7f
312 axis = (word0 >> 23) & 0x01
313 lane = (word0 >> 24) & 0x1f # 1..2 for scints, 8..20 for RPCs (=readout-board slot - 7)
314 flag = (word0 >> 30) & 0x03 # 1 for RPCs, 2 for scints
315 tdc = (word1 >> 16) & 0x07ff
316 isRPC = (flag == 1)
317 electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
318 if electId in self.electIdToModuleId: # BKLM mapped-channel histograms
319 if isRPC:
320 tCal = int(tdc - trigCtime) & 0x03ff # in ns, range is 0..1023
321 self.hist_mappedRPCTimeCal.Fill(tCal)
322 self.hist_mappedRPCTimeCalByWindow.Fill(tCal)
323
324 # for normalization of the hit-counter dictionaries, now that we know that this is a valid event
325
326 self.dict_nEvents[self.windowValue] += 1
327
328 # Process the KLMHit2ds
329
330 self.dict_nHit2ds[self.windowValue] += len(hit2ds)
331 for hit2d in hit2ds:
332 key = hit2d.getModuleID()
333 fb = (key & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
334 x = hit2d.getGlobalPositionX()
335 y = hit2d.getGlobalPositionY()
336 if fb == 0: # backward
337 self.hist_occupancyBackwardXY.Fill(x, y)
338 else: # forward
339 self.hist_occupancyForwardXY.Fill(x, y)
340 self.hist_occupancyXYByWindow.Fill(x, y)
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72

◆ initialize()

def initialize (   self)
Handle job initialization: fill the mapping database, create histograms

Definition at line 93 of file EventInspectorLookback.py.

93 def initialize(self):
94 """Handle job initialization: fill the mapping database, create histograms"""
95
96 expRun = f'e{int(self.exp):02d}r{int(self.run)}: '
97
98
99 self.electIdToModuleId = bklmDB.fillDB()
100
101 self.sectorFBToDC = [11, 15, 2, 6, 10, 14, 3, 7, 9, 13, 0, 4, 8, 12, 1, 5]
102
103 self.dcToSectorFB = [10, 14, 2, 6, 11, 15, 3, 7, 12, 8, 4, 0, 13, 9, 5, 1]
104
105 self.histogramFile = ROOT.TFile.Open(self.histName, "RECREATE")
106 # All histograms/scatterplots in the output file will show '# of events' only
107 ROOT.gStyle.SetOptStat(10)
108
109 # create the rawKLM histograms
110
111
112 self.hist_mappedRPCTimeCal = ROOT.TH1F(
113 'mappedRPCTimeCal', expRun + 'RPC time distribution;t - t(trigger) (ns)', 256, -0.5, 1023.5)
114
115 self.dict_mappedRPCTimeCalByWindow = {}
116
117 self.dict_nRawKLMs = {}
118 for window in range(self.windowMinValue, self.windowMaxValue+1, self.windowStepValue):
119 label = f"mappedRPCTimeCalByWindow{window:04x}"
120 title = f"{expRun}RPC time distribution for lookback-window {self.windowMode} = {window:04x}" + \
121 ";t - t(trigger) (ns)"
122 self.dict_mappedRPCTimeCalByWindow[window] = ROOT.TH1F(label, title, 256, -0.5, 1023.5)
123 self.dict_nRawKLMs[window] = 0
124
125 self.hist_mappedRPCTimeCalByWindow = self.dict_mappedRPCTimeCalByWindow[self.windowMinValue]
126
127 # Create the KLMHit2d-related histograms
128
129
130 self.hist_occupancyForwardXY = ROOT.TH2F('occupancyForwardXY',
131 expRun + 'Forward xy occupancy;x(cm);y(cm)',
132 230, -345.0, 345.0, 230, -345.0, 345.0)
133
134 self.hist_occupancyBackwardXY = ROOT.TH2F('occupancyBackwardXY',
135 expRun + 'Backward xy occupancy;x(cm);y(cm)',
136 230, -345.0, 345.0, 230, -345.0, 345.0)
137
138 self.hist_occupancyXYByWindow = self.dict_occupancyXYByWindow[self.windowMinValue]
139
140
141 self.dict_occupancyXYByWindow = {}
142
143 self.dict_nHit2ds = {}
144 for window in range(self.windowMinValue, self.windowMaxValue+1, self.windowStepValue):
145 label = f"occupancyXYByWindow{window:04x}"
146 title = f"{expRun}Forward xy occupancy for lookback-window {self.windowMode} = {window:04x};x(cm);y(cm)"
147 self.dict_occupancyXYByWindow[window] = ROOT.TH2F(label, title, 230, -345.0, 345.0, 230, -345.0, 345.0)
148 self.dict_nHit2ds[window] = 0
149
150
151 self.dict_nEvents = {}
152 for window in range(self.windowMinValue, self.windowMaxValue+1, self.windowStepValue):
153 self.dict_nEvents[window] = 0
154
155
156 self.windowValue = -1
157
def fillDB()
Definition: bklmDB.py:16

◆ terminate()

def terminate (   self)
Handle job termination: draw histograms, close output files

Definition at line 158 of file EventInspectorLookback.py.

158 def terminate(self):
159 """Handle job termination: draw histograms, close output files"""
160
161 canvas = ROOT.TCanvas("canvas", self.pdfName, 1600, 1600)
162 title = f'{self.pdfName}['
163 canvas.SaveAs(title)
164 canvas.Clear()
165 canvas.Divide(1, 1)
166 canvas.GetPad(0).SetGrid(1, 1)
167 canvas.GetPad(0).Update()
168 nWindows = int((self.windowMaxValue - self.windowMinValue) / self.windowStepValue) + 1
169 hist_nEvents = ROOT.TH1F('nEvents',
170 f'Number of events;Lookback-window {self.windowMode} value',
171 nWindows, self.windowMinValue, self.windowMaxValue+self.windowStepValue)
172 hist_nEvents.SetStats(False)
173 hist_nEvents.SetMinimum(0)
174 values = array.array('d', [0]) # dummy value for the histogram's underflow bin
175 for key in self.dict_nEvents:
176 values.append(self.dict_nEvents[key])
177 values.append(0) # dummy value for the histogram's overflow bin
178 hist_nEvents.SetContent(values)
179 hist_nEvents.Draw("HIST")
180 canvas.Print(self.pdfName, f"Title:{hist_nEvents.GetName()}")
181 hist_nRawKLMs = ROOT.TH1F('nRawKLMs',
182 f'Mean number of RawKLM hits per event;Lookback-window {self.windowMode} value',
183 nWindows, self.windowMinValue, self.windowMaxValue+self.windowStepValue)
184 hist_nRawKLMs.SetStats(False)
185 hist_nRawKLMs.SetMinimum(0)
186 ratios = array.array('d', [0]) # dummy ratio for the histogram's underflow bin
187 errors = array.array('d', [0]) # dummy error for the histogram's underflow bin
188 for key in self.dict_nRawKLMs:
189 numerator = self.dict_nRawKLMs[key]
190 denominator = float(self.dict_nEvents[key])
191 if denominator > 0:
192 ratio = numerator / denominator
193 ratios.append(ratio)
194 errors.append(math.sqrt(ratio * (ratio + 1.0) / denominator)) # avoid 1/numerator
195 else:
196 ratios.append(0)
197 errors.append(0)
198 ratios.append(0) # dummy ratio for the histogram's overflow bin
199 errors.append(0) # dummy error for the histogram's overflow bin
200 hist_nRawKLMs.SetContent(ratios)
201 hist_nRawKLMs.SetError(errors)
202 hist_nRawKLMs.Draw("E0 X0 L")
203 hist_nRawKLMs.Draw("HIST SAME")
204 canvas.Print(self.pdfName, f"Title:{hist_nRawKLMs.GetName()}")
205 hist_nHit2ds = ROOT.TH1F('nKLMHit2ds',
206 f'Mean number of KLMHit2ds per event;Lookback-window {self.windowMode} value',
207 nWindows, self.windowMinValue, self.windowMaxValue+self.windowStepValue)
208 hist_nHit2ds.SetStats(False)
209 hist_nHit2ds.SetMinimum(0)
210 ratios = array.array('d', [0]) # dummy ratio for the histogram's underflow bin
211 errors = array.array('d', [0]) # dummy error for the histogram's underflow bin
212 for key in self.dict_nHit2ds:
213 numerator = self.dict_nHit2ds[key]
214 denominator = float(self.dict_nEvents[key])
215 if denominator > 0:
216 ratio = numerator / denominator
217 ratios.append(ratio)
218 errors.append(math.sqrt(ratio * (ratio + 1.0) / denominator)) # avoid 1/numerator
219 else:
220 ratios.append(0)
221 errors.append(0)
222 ratios.append(0) # dummy ratio for the histogram's overflow bin
223 errors.append(0) # dummy error for the histogram's overflow bin
224 hist_nHit2ds.SetContent(ratios)
225 hist_nHit2ds.SetError(errors)
226 hist_nHit2ds.Draw("E0 X0 L")
227 hist_nHit2ds.Draw("HIST SAME")
228 canvas.Print(self.pdfName, f"Title:{hist_nHit2ds.GetName()}")
229
230 self.hist_mappedRPCTimeCal.Draw()
231 canvas.Print(self.pdfName, f"Title:{self.hist_mappedRPCTimeCal.GetName()}")
232 for key in self.dict_mappedRPCTimeCalByWindow:
233 theHist = self.dict_mappedRPCTimeCalByWindow[key]
234 theHist.Draw()
235 canvas.Print(self.pdfName, f"Title:{theHist.GetName()}")
236
237 self.hist_occupancyBackwardXY.Draw("colz")
238 canvas.Print(self.pdfName, f"Title:{self.hist_occupancyBackwardXY.GetName()}")
239 self.hist_occupancyForwardXY.Draw("colz")
240 canvas.Print(self.pdfName, f"Title:{self.hist_occupancyForwardXY.GetName()}")
241 for key in self.dict_occupancyXYByWindow:
242 theHist = self.dict_occupancyXYByWindow[key]
243 theHist.Draw("colz")
244 lastTitle = f"Title:{theHist.GetName()}"
245 canvas.Print(self.pdfName, lastTitle)
246 pdfNameLast = f'{self.pdfName}]'
247 canvas.Print(pdfNameLast, lastTitle)
248 self.histogramFile.Write()
249 self.histogramFile.Close()
250 print('Goodbye')
251

Member Data Documentation

◆ BKLM_ID

int BKLM_ID = 0x07000000
static

COPPER base identifier for BKLM readout.

Definition at line 27 of file EventInspectorLookback.py.

◆ BKLM_LAYER_BIT

int BKLM_LAYER_BIT = 7
static

bit position for layer-1 [0..14]; 0 is innermost

Definition at line 35 of file EventInspectorLookback.py.

◆ BKLM_LAYER_MASK

tuple BKLM_LAYER_MASK = (15 << BKLM_LAYER_BIT)
static

bit mask for layer-1 [0..15]; 0 is innermost and 14 is outermost

Definition at line 47 of file EventInspectorLookback.py.

◆ BKLM_MAXSTRIP_BIT

int BKLM_MAXSTRIP_BIT = 15
static

bit position for maxStrip-1 [0..47]

Definition at line 41 of file EventInspectorLookback.py.

◆ BKLM_MAXSTRIP_MASK

tuple BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
static

bit mask for maxStrip-1 [0..47]

Definition at line 53 of file EventInspectorLookback.py.

◆ BKLM_MODULEID_MASK

tuple BKLM_MODULEID_MASK = (BKLM_SECTION_MASK | BKLM_SECTOR_MASK | BKLM_LAYER_MASK)
static

bit mask for unique module identifier (end, sector, layer)

Definition at line 55 of file EventInspectorLookback.py.

◆ BKLM_PLANE_BIT

int BKLM_PLANE_BIT = 6
static

bit position for plane-1 [0..1]; 0 is inner-plane

Definition at line 33 of file EventInspectorLookback.py.

◆ BKLM_PLANE_MASK

tuple BKLM_PLANE_MASK = (1 << BKLM_PLANE_BIT)
static

bit mask for plane-1 [0..1]; 0 is inner-plane

Definition at line 45 of file EventInspectorLookback.py.

◆ BKLM_SECTION_BIT

int BKLM_SECTION_BIT = 14
static

bit position for section [0..1]; forward is 0

Definition at line 39 of file EventInspectorLookback.py.

◆ BKLM_SECTION_MASK

tuple BKLM_SECTION_MASK = (1 << BKLM_SECTION_BIT)
static

bit mask for section [0..1]; forward is 0

Definition at line 51 of file EventInspectorLookback.py.

◆ BKLM_SECTOR_BIT

int BKLM_SECTOR_BIT = 11
static

bit position for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis

Definition at line 37 of file EventInspectorLookback.py.

◆ BKLM_SECTOR_MASK

tuple BKLM_SECTOR_MASK = (7 << BKLM_SECTOR_BIT)
static

bit mask for sector-1 [0..7]; 0 is on the +x axis and 2 is on the +y axis

Definition at line 49 of file EventInspectorLookback.py.

◆ BKLM_STRIP_BIT

int BKLM_STRIP_BIT = 0
static

bit position for strip-1 [0..47]

Definition at line 31 of file EventInspectorLookback.py.

◆ BKLM_STRIP_MASK

int BKLM_STRIP_MASK = 0x3f
static

bit mask for strip-1 [0..47]

Definition at line 43 of file EventInspectorLookback.py.

◆ dcToSectorFB

dcToSectorFB

map for data concentrator -> sectorFB

Definition at line 103 of file EventInspectorLookback.py.

◆ dict_mappedRPCTimeCalByWindow

dict_mappedRPCTimeCalByWindow

dictionary of histograms of RPC TDC - trigger value, keyed by lookback-window value

Definition at line 115 of file EventInspectorLookback.py.

◆ dict_nEvents

dict_nEvents

dictionary of the number of events for each lookback-window value, for normalization

Definition at line 151 of file EventInspectorLookback.py.

◆ dict_nHit2ds

dict_nHit2ds

dictionary of the number of KLMHit2ds for each lookback-window value

Definition at line 143 of file EventInspectorLookback.py.

◆ dict_nRawKLMs

dict_nRawKLMs

dictionary of the number of RawKLM hits for each lookback-window value

Definition at line 117 of file EventInspectorLookback.py.

◆ dict_occupancyXYByWindow

dict_occupancyXYByWindow

dictionary of scatterplots of end view of forward BKLM, keyed by lookback-window value

Definition at line 141 of file EventInspectorLookback.py.

◆ EKLM_ID

int EKLM_ID = 0x08000000
static

COPPER base identifier for EKLM readout.

Definition at line 29 of file EventInspectorLookback.py.

◆ electIdToModuleId

electIdToModuleId

readout <-> detector map (from the information retrieved from the conditions database)

Definition at line 99 of file EventInspectorLookback.py.

◆ exp

exp

internal copy of experiment number

Definition at line 74 of file EventInspectorLookback.py.

◆ hist_mappedRPCTimeCal

hist_mappedRPCTimeCal

histogram of RPC TDC - trigger value

Definition at line 112 of file EventInspectorLookback.py.

◆ hist_mappedRPCTimeCalByWindow

hist_mappedRPCTimeCalByWindow

reference to the RPC-time histogram for the currevent value of the lookback window parameter

Definition at line 125 of file EventInspectorLookback.py.

◆ hist_occupancyBackwardXY

hist_occupancyBackwardXY

scatterplot of end view of backward BKLM for all KLMHit2ds

Definition at line 134 of file EventInspectorLookback.py.

◆ hist_occupancyForwardXY

hist_occupancyForwardXY

scatterplot of end view of forward BKLM for all KLMHit2ds

Definition at line 130 of file EventInspectorLookback.py.

◆ hist_occupancyXYByWindow

hist_occupancyXYByWindow

reference to the xy scatterplot for the currevent value of the lookback window parameter

Definition at line 138 of file EventInspectorLookback.py.

◆ histName

histName

internal copy of the pathname of the output histogram ROOT file

Definition at line 78 of file EventInspectorLookback.py.

◆ histogramFile

histogramFile

Output ROOT TFile that will contain the histograms/scatterplots.

Definition at line 105 of file EventInspectorLookback.py.

◆ pdfName

pdfName

internal copy of the pathname of the output histogram PDF file

Definition at line 80 of file EventInspectorLookback.py.

◆ run

run

internal copy of run number

Definition at line 76 of file EventInspectorLookback.py.

◆ sectorFBToDC

sectorFBToDC

map for sectorFB -> data concentrator

Definition at line 101 of file EventInspectorLookback.py.

◆ windowMaxValue

windowMaxValue

highest observed lookback-window value

Definition at line 88 of file EventInspectorLookback.py.

◆ windowMinValue

windowMinValue

highest observed lookback-window value

Definition at line 86 of file EventInspectorLookback.py.

◆ windowMode

windowMode

window mode as a string for histogram labels/titles

Definition at line 84 of file EventInspectorLookback.py.

◆ windowStepValue

windowStepValue

lookback-window value step

Definition at line 90 of file EventInspectorLookback.py.

◆ windowValue

windowValue

cached value of the lookback-window value, to avoid unnecessary reassignments-to-same-value in event()

Definition at line 156 of file EventInspectorLookback.py.


The documentation for this class was generated from the following file: