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

Public Member Functions

def __init__ (self, exp, run, eventPdfName, maxDisplays, minRPCHits, minMuidHits)
 
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
 
 eventPdfName
 internal copy of the pathname of the output event-display PDF file
 
 maxDisplays
 internal copy of the maximum number of event displays to write
 
 minRPCHits
 internal copy of the minimum number of RPC KLMHit2ds in any sector for event display
 
 minMuidHits
 internal copy of the minimum number of MuidHits in the event for event display
 
 eventCounter
 event counter (needed for PDF table of contents' ordinal event#)
 
 eventDisplays
 event-display counter
 
 lastTitle
 title of the last-drawn event display (needed for PDF table of contents' last event)
 
 eventCanvas
 TCanvas on which event displays will be drawn.
 
 cosine
 table of cosines for the BKLM-sector normals
 
 sine
 table of sines for the BKLM-sector normals
 
 hist_XY
 blank scatterplot to define the bounds of the BKLM end view
 
 hist_ZY
 blank scatterplot to define the bounds of the BKLM side view
 
 hist_XYS
 list of blank scatterplots to define the per-sector bounds of the BKLM end view
 
 hist_ZYS
 blank scatterplot to define the per-sector bounds of the rotated BKLM side view
 
 bklmXY
 list of line-segment (x,y) points for the BKLM end view
 
 bklmZY
 list of line-segment (z,y) points for the BKLM side view
 
 bklmZYL
 list of line-segment (z,y) points for the BKLM sector's zoomed and rotated side view
 
 electIdToModuleId
 readout <-> detector map (from the information retrieved from the conditions database)
 
 sectorFBToDC
 map for sectorFB -> data concentrator
 
 dcToSectorFB
 map for data concentrator -> sectorFB
 
 t0Cal
 RPC-time calibration adjustment (ns) for rawKLMs.
 
 t0Cal2d
 RPC-time calibration adjustment (ns) for KLMHit2ds.
 
 ct0Cal
 scint-ctime calibration adjustment (ns) for rawKLMs
 
 ct0Cal2d
 scint-ctime calibration adjustment (ns) for KLMHit2ds
 
 t0RPC
 per-sector variations in RPC-time calibration adjustment (ns) for rawKLMs
 
 ct0Scint
 per-sector variations in scint-ctime calibration adjustment (ns) for rawKLMs
 

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

Draw BKLM event displays from KLMHit2ds, ExtHits, and MuidHits.

Definition at line 28 of file EventDisplayer.py.

Constructor & Destructor Documentation

◆ __init__()

def __init__ (   self,
  exp,
  run,
  eventPdfName,
  maxDisplays,
  minRPCHits,
  minMuidHits 
)
Constructor

Arguments:
    exp (str): formatted experiment number
    run (str): formatter run number
    eventPdfName (str): path name of the output event-display PDF file
    maxDisplays (int): max # of events displays to write
    minRPCHits (int): min # of RPC KLMHit2ds in any sector for event display
    minMuidHits (int): min # of MuidHits in the event for event display

Definition at line 62 of file EventDisplayer.py.

62 def __init__(self, exp, run, eventPdfName, maxDisplays, minRPCHits, minMuidHits):
63 """Constructor
64
65 Arguments:
66 exp (str): formatted experiment number
67 run (str): formatter run number
68 eventPdfName (str): path name of the output event-display PDF file
69 maxDisplays (int): max # of events displays to write
70 minRPCHits (int): min # of RPC KLMHit2ds in any sector for event display
71 minMuidHits (int): min # of MuidHits in the event for event display
72 """
73 super().__init__()
74
75 self.exp = exp
76
77 self.run = run
78
79 self.eventPdfName = eventPdfName
80
81 self.maxDisplays = maxDisplays
82
83 self.minRPCHits = minRPCHits
84
85 self.minMuidHits = minMuidHits
86
87 self.eventCounter = 0
88
89 self.eventDisplays = 0
90
91 self.lastTitle = ''
92

Member Function Documentation

◆ beginRun()

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

Definition at line 320 of file EventDisplayer.py.

320 def beginRun(self):
321 """Handle begin of run: print diagnostic message"""
322 EventMetaData = Belle2.PyStoreObj('EventMetaData')
323 print('beginRun', EventMetaData.getRun())
324
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67

◆ endRun()

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

Definition at line 325 of file EventDisplayer.py.

325 def endRun(self):
326 """Handle end of run: print diagnostic message"""
327 EventMetaData = Belle2.PyStoreObj('EventMetaData')
328 print('endRun', EventMetaData.getRun())
329

◆ event()

def event (   self)
Process one event: (optionally) draw event display

Definition at line 330 of file EventDisplayer.py.

330 def event(self):
331 """Process one event: (optionally) draw event display"""
332
333 super().return_value(self.eventDisplays < self.maxDisplays)
334
335 if self.eventDisplays >= self.maxDisplays:
336 return
337
338 self.eventCounter += 1
339 EventMetaData = Belle2.PyStoreObj('EventMetaData')
340 event = EventMetaData.getEvent()
341 rawklms = Belle2.PyStoreArray('RawKLMs') # to determine if KLMHit2d is prompt
342 hit2ds = Belle2.PyStoreArray('KLMHit2ds')
343 exthits = Belle2.PyStoreArray('ExtHits')
344 muidhits = Belle2.PyStoreArray('MuidHits')
345
346 rawFb = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
347 rawSector = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
348 rawLayer = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
349 rawPlane = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
350 rawStrip = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
351 rawCtime = [[], [], [], [], [], [], [], [], [], [], [], [], [], [], [], []]
352
353 for copper in range(0, len(rawklms)):
354 rawklm = rawklms[copper]
355 nodeID = rawklm.GetNodeID(0) - self.BKLM_ID
356 if nodeID >= self.EKLM_ID - self.BKLM_ID:
357 nodeID = nodeID - (self.EKLM_ID - self.BKLM_ID) + 4
358 if (nodeID < 0) or (nodeID > 4): # skip EKLM nodes
359 continue
360 for finesse in range(0, 4):
361 dc = (finesse << 2) + copper
362 n = rawklm.GetDetectorNwords(0, finesse)
363 bufSlot = rawklm.GetDetectorBuffer(0, finesse)
364 if n > 0:
365 n = n >> 1
366 for j in range(0, n):
367 word0 = bufSlot[j * 2]
368 word1 = bufSlot[j * 2 + 1] # noqa (F841) kept for completeness
369 ctime = word0 & 0xffff
370 channel = (word0 >> 16) & 0x7f
371 axis = (word0 >> 23) & 0x01
372 lane = (word0 >> 24) & 0x1f # crate's slot number
373 flag = (word0 >> 30) & 0x03 # noqa (F841) kept for completeness
374 fb = -1
375 sector = -1
376 layer = -1
377 plane = -1
378 strip = -1
379 electId = (channel << 12) | (axis << 11) | (lane << 6) | (finesse << 4) | nodeID
380 if electId in self.electIdToModuleId:
381 moduleId = self.electIdToModuleId[electId]
382 fb = (moduleId & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
383 sector = (moduleId & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
384 layer = (moduleId & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
385 plane = (moduleId & self.BKLM_PLANE_MASK) >> self.BKLM_PLANE_BIT
386 strip = (moduleId & self.BKLM_STRIP_MASK) >> self.BKLM_STRIP_BIT
387 rawFb[dc].append(fb)
388 rawSector[dc].append(sector)
389 rawLayer[dc].append(layer)
390 rawPlane[dc].append(plane)
391 rawStrip[dc].append(strip)
392 rawCtime[dc].append(ctime)
393
394 tCal2d = []
395 for hit2d in hit2ds:
396 tCal2d.append(hit2d.getTime())
397
398 # Process the ExtHits for event display
399
400 xold = 0
401 yold = 0
402 zold = 0
403 sold = 0
404 epositions = []
405 sumx = 0
406 sumy = 0
407 sumz = 0
408 sumn = 0
409 for exthit in exthits:
410 sector = (exthit.getCopyID() & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
411 extPosition = exthit.getPosition()
412 x = extPosition[0]
413 y = extPosition[1]
414 z = extPosition[2]
415 dx = x - xold
416 dy = y - yold
417 dz = z - zold
418 if (dx * dx + dy * dy + dz * dz > 36) and (sumn > 0):
419 eposition = [sumx / sumn, sumy / sumn, sumz / sumn, sold]
420 epositions.append(eposition)
421 sumx = 0
422 sumy = 0
423 sumz = 0
424 sumn = 0
425 else:
426 if sumn == 0:
427 xold = x
428 yold = y
429 zold = z
430 sold = sector
431 sumx = sumx + x
432 sumy = sumy + y
433 sumz = sumz + z
434 sumn = sumn + 1
435 if sumn > 0:
436 eposition = [sumx / sumn, sumy / sumn, sumz / sumn, sold]
437 epositions.append(eposition)
438 extXYGraph = ROOT.TGraph()
439 extXYGraph.SetMarkerColor(30)
440 extXYGraph.SetMarkerSize(2.25)
441 extXYGraph.SetMarkerStyle(21)
442 extZYGraph = ROOT.TGraph()
443 extZYGraph.SetMarkerColor(30)
444 extZYGraph.SetMarkerSize(2.25)
445 extZYGraph.SetMarkerStyle(21)
446 extZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
447 for sector in range(0, 8):
448 extZYSGraph[sector] = ROOT.TGraph()
449 extZYSGraph[sector].SetMarkerColor(30)
450 extZYSGraph[sector].SetMarkerSize(2.25)
451 extZYSGraph[sector].SetMarkerStyle(21)
452 j = -1
453 for eposition in epositions:
454 j = j + 1
455 x = eposition[0]
456 y = eposition[1]
457 z = eposition[2]
458 extXYGraph.SetPoint(j, x, y)
459 extZYGraph.SetPoint(j, z, y)
460 sector = int(eposition[3])
461 nPoint = extZYSGraph[sector].GetN()
462 extZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosine[sector] + y * self.sine[sector]))
463
464 # Process the MuidHits for event display
465
466 zMuids = [0, 0, 0, 0, 0, 0, 0, 0]
467 nMuids = [0, 0, 0, 0, 0, 0, 0, 0]
468 muidXYGraph = ROOT.TGraph()
469 muidXYGraph.SetMarkerColor(5)
470 muidXYGraph.SetMarkerSize(2.0)
471 muidXYGraph.SetMarkerStyle(20)
472 muidZYGraph = ROOT.TGraph()
473 muidZYGraph.SetMarkerColor(5)
474 muidZYGraph.SetMarkerSize(2.0)
475 muidZYGraph.SetMarkerStyle(20)
476 muidZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
477 for sector in range(0, 8):
478 muidZYSGraph[sector] = ROOT.TGraph()
479 muidZYSGraph[sector].SetMarkerColor(5)
480 muidZYSGraph[sector].SetMarkerSize(2.0)
481 muidZYSGraph[sector].SetMarkerStyle(20)
482 j = -1
483 for muidhit in muidhits:
484 j += 1
485 muidPosition = muidhit.getExtPosition()
486 x = muidPosition[0]
487 y = muidPosition[1]
488 z = muidPosition[2]
489 muidXYGraph.SetPoint(j, x, y)
490 muidZYGraph.SetPoint(j, z, y)
491 sector = muidhit.getSector()
492 nPoint = muidZYSGraph[sector].GetN()
493 muidZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosine[sector] + y * self.sine[sector]))
494 if nMuids[sector] == 0:
495 zMuids[sector] = z
496 nMuids[sector] += 1
497
498 # Process the KLMHit2ds for event display
499
500 rpcHitCount = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
501 promptXYGraph = ROOT.TGraph()
502 promptXYGraph.SetMarkerColor(4)
503 promptXYGraph.SetMarkerSize(2.0)
504 promptXYGraph.SetMarkerStyle(29)
505 promptZYGraph = ROOT.TGraph()
506 promptZYGraph.SetMarkerColor(4)
507 promptZYGraph.SetMarkerSize(2.0)
508 promptZYGraph.SetMarkerStyle(29)
509 bkgdXYGraph = ROOT.TGraph()
510 bkgdXYGraph.SetMarkerColor(2)
511 bkgdXYGraph.SetMarkerSize(2.0)
512 bkgdXYGraph.SetMarkerStyle(29)
513 bkgdZYGraph = ROOT.TGraph()
514 bkgdZYGraph.SetMarkerColor(2)
515 bkgdZYGraph.SetMarkerSize(2.0)
516 bkgdZYGraph.SetMarkerStyle(29)
517 promptZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
518 bkgdZYSGraph = [0, 0, 0, 0, 0, 0, 0, 0]
519 for sector in range(0, 8):
520 promptZYSGraph[sector] = ROOT.TGraph()
521 promptZYSGraph[sector].SetMarkerColor(4)
522 promptZYSGraph[sector].SetMarkerSize(2.0)
523 promptZYSGraph[sector].SetMarkerStyle(29)
524 bkgdZYSGraph[sector] = ROOT.TGraph()
525 bkgdZYSGraph[sector].SetMarkerColor(2)
526 bkgdZYSGraph[sector].SetMarkerSize(2.0)
527 bkgdZYSGraph[sector].SetMarkerStyle(29)
528 jPrompt = -1
529 jBkgd = -1
530 for hit2d in hit2ds:
531 key = hit2d.getModuleID()
532 layer = (key & self.BKLM_LAYER_MASK) >> self.BKLM_LAYER_BIT
533 sector = (key & self.BKLM_SECTOR_MASK) >> self.BKLM_SECTOR_BIT
534 fb = (key & self.BKLM_SECTION_MASK) >> self.BKLM_SECTION_BIT
535 phiStripMin = hit2d.getPhiStripMin() - 1
536 phiStripMax = hit2d.getPhiStripMax() - 1
537 zStripMin = hit2d.getZStripMin() - 1
538 zStripMax = hit2d.getZStripMax() - 1
539 sectorFB = sector if fb == 0 else sector + 8
540 if layer >= 2:
541 rpcHitCount[sectorFB] += 1
542 dc = self.sectorFBToDC[sectorFB]
543 copper = dc & 0x03
544 finesse = dc >> 2
545 n = rawklms[copper].GetDetectorNwords(0, finesse) >> 1
546 trigCtime = (rawklms[copper].GetTTCtime(0) & 0x07ffffff) << 3
547 ctDiffMax = 99999
548 tCal = -1
549 jZ = -1
550 jPhi = -1
551 ctZ = 0
552 ctPhi = 0
553 for j in range(0, n):
554 if layer != rawLayer[dc][j]:
555 continue
556 if sector != rawSector[dc][j]:
557 continue
558 if fb != rawFb[dc][j]:
559 continue
560 strip = rawStrip[dc][j]
561 plane = rawPlane[dc][j]
562 if plane == 0: # it's a z strip
563 if strip < zStripMin:
564 continue
565 if strip > zStripMax:
566 continue
567 ctZ = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
568 jZ = j
569 else: # it's a phi strip
570 if strip < phiStripMin:
571 continue
572 if strip > phiStripMax:
573 continue
574 ctPhi = rawCtime[dc][j] << 3 # in ns, range is only 8 bits in SCROD (??)
575 jPhi = j
576 if (jZ >= 0) and (jPhi >= 0):
577 if layer < 2: # it's a scint layer
578 if abs(ctZ - ctPhi) > 40:
579 continue
580 ct = int((ctZ + ctPhi) * 0.5 - trigCtime - self.ct0Scint[sectorFB]) & 0x3ff
581 if abs(ct - self.ct0Cal) < ctDiffMax:
582 ctDiffMax = int(abs(ct - self.ct0Cal))
583 tCal = ct
584 if ctDiffMax == 0:
585 break
586 else: # it's an RPC layer
587 tCal = ((int(hit2d.getTime()) - trigCtime) & 0x03ff) - self.t0RPC[sectorFB] - 0.75 * jPhi - 0.75 * jZ
588 break
589 x = hit2d.getGlobalPositionX()
590 y = hit2d.getGlobalPositionY()
591 z = hit2d.getGlobalPositionZ()
592 isPromptHit = False
593 if layer < 2:
594 if abs(tCal - self.ct0Cal2d) < 20:
595 isPromptHit = True
596 else:
597 if abs(tCal - self.t0Cal2d) < 20:
598 isPromptHit = True
599 if isPromptHit:
600 jPrompt += 1
601 promptXYGraph.SetPoint(jPrompt, x, y)
602 promptZYGraph.SetPoint(jPrompt, z, y)
603 nPoint = promptZYSGraph[sector].GetN()
604 promptZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosine[sector] + y * self.sine[sector]))
605 else:
606 jBkgd += 1
607 bkgdXYGraph.SetPoint(jBkgd, x, y)
608 bkgdZYGraph.SetPoint(jBkgd, z, y)
609 nPoint = bkgdZYSGraph[sector].GetN()
610 bkgdZYSGraph[sector].SetPoint(nPoint, z, abs(x * self.cosine[sector] + y * self.sine[sector]))
611
612 hasEnoughRPCHits = False
613 for count in rpcHitCount:
614 if count > self.minRPCHits:
615 hasEnoughRPCHits = True
616 break
617 if hasEnoughRPCHits and (len(muidhits) > self.minMuidHits):
618 self.eventDisplays += 1
619 title = f'e{int(self.exp):02d}r{int(self.run)}: event {event}'
620 self.hist_XY.SetTitle(title)
621 self.hist_ZY.SetTitle(title)
622 self.eventCanvas.cd(1)
623 self.hist_XY.Draw()
624 for g in self.bklmXY:
625 g.Draw("L")
626 if extXYGraph.GetN() > 0:
627 extXYGraph.Draw("P")
628 if muidXYGraph.GetN() > 0:
629 muidXYGraph.Draw("P")
630 if bkgdXYGraph.GetN() > 0:
631 bkgdXYGraph.Draw("P")
632 if promptXYGraph.GetN() > 0:
633 promptXYGraph.Draw("P")
634 self.eventCanvas.cd(2)
635 self.hist_ZY.Draw()
636 for g in self.bklmZY:
637 g.Draw("L")
638 if extZYGraph.GetN() > 0:
639 extZYGraph.Draw("P")
640 if muidZYGraph.GetN() > 0:
641 muidZYGraph.Draw("P")
642 if bkgdZYGraph.GetN() > 0:
643 bkgdZYGraph.Draw("P")
644 if promptZYGraph.GetN() > 0:
645 promptZYGraph.Draw("P")
646 self.lastTitle = f"Title:E{event} (#{self.eventCounter})"
647 self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
648 for sector in range(0, 8):
649 if nMuids[sector] > 0:
650 title = f'e{int(self.exp):02d}r{int(self.run)}: event {event} sector {sector}'
651 self.hist_XYS[sector].SetTitle(title)
652 self.eventCanvas.cd(1)
653 self.hist_XYS[sector].Draw()
654 for g in self.bklmXY:
655 g.Draw("L")
656 if extXYGraph.GetN() > 0:
657 extXYGraph.Draw("P")
658 if muidXYGraph.GetN() > 0:
659 muidXYGraph.Draw("P")
660 if bkgdXYGraph.GetN() > 0:
661 bkgdXYGraph.Draw("P")
662 if promptXYGraph.GetN() > 0:
663 promptXYGraph.Draw("P")
664 self.eventCanvas.cd(2)
665 z0 = zMuids[sector]
666 self.hist_ZYS.SetTitle(title)
667 self.hist_ZYS.SetBins(10, z0 - 150.0, z0 + 150.0, 10, 125.0, 425.0)
668 self.hist_ZYS.Draw()
669 for g in self.bklmZY:
670 g.Draw("L")
671 for g in self.bklmZYL:
672 g.Draw("L")
673 if extZYSGraph[sector].GetN() > 0:
674 extZYSGraph[sector].Draw("P")
675 if muidZYSGraph[sector].GetN() > 0:
676 muidZYSGraph[sector].Draw("P")
677 if bkgdZYSGraph[sector].GetN() > 0:
678 bkgdZYSGraph[sector].Draw("P")
679 if promptZYSGraph[sector].GetN() > 0:
680 promptZYSGraph[sector].Draw("P")
681 self.lastTitle = f"Title:E{event} sector {sector}"
682 self.eventCanvas.Print(self.eventPdfName, self.lastTitle)
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72

◆ initialize()

def initialize (   self)
Handle job initialization: fill the mapping database, create histograms, open the event-display file

Definition at line 93 of file EventDisplayer.py.

93 def initialize(self):
94 """Handle job initialization: fill the mapping database, create histograms, open the event-display file"""
95
96 # Open the output PDF file for event displays
97
98
99 self.eventCanvas = ROOT.TCanvas("eventCanvas", self.eventPdfName, 3200, 1600)
100 title = f'{self.eventPdfName}['
101 self.eventCanvas.SaveAs(title)
102 self.eventCanvas.Clear()
103 self.eventCanvas.Divide(2, 1)
104
105 # Create the boilerplate for the end- and side-views of the event display
106
107
108 self.cosine = [0, 0, 0, 0, 0, 0, 0, 0]
109
110 self.sine = [0, 0, 0, 0, 0, 0, 0, 0]
111 for sector in range(0, 8):
112 phi = math.pi * sector / 4
113 self.cosine[sector] = math.cos(phi)
114 self.sine[sector] = math.sin(phi)
115
116 self.hist_XY = ROOT.TH2F('XY', ' ;x;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
117 self.hist_XY.SetStats(False)
118
119 self.hist_ZY = ROOT.TH2F('ZY', ' ;z;y', 10, -345.0, 345.0, 10, -345.0, 345.0)
120 self.hist_ZY.SetStats(False)
121 # 300x300 cm^2 grid for each octant
122 u1 = 65
123 u2 = 365
124 u3 = 100
125 u4 = 400
126 u5 = 150
127
128 self.hist_XYS = [0, 0, 0, 0, 0, 0, 0, 0]
129 self.hist_XYS[0] = ROOT.TH2F('XYS0', ' ;x;y', 10, +u3, +u4, 10, -u5, +u5)
130 self.hist_XYS[0].SetStats(False)
131 self.hist_XYS[1] = ROOT.TH2F('XYS1', ' ;x;y', 10, +u1, +u2, 10, +u1, +u2)
132 self.hist_XYS[1].SetStats(False)
133 self.hist_XYS[2] = ROOT.TH2F('XYS2', ' ;x;y', 10, -u5, +u5, 10, +u3, +u4)
134 self.hist_XYS[2].SetStats(False)
135 self.hist_XYS[3] = ROOT.TH2F('XYS3', ' ;x;y', 10, -u2, -u1, 10, +u1, +u2)
136 self.hist_XYS[3].SetStats(False)
137 self.hist_XYS[4] = ROOT.TH2F('XYS4', ' ;x;y', 10, -u4, -u3, 10, -u5, +u5)
138 self.hist_XYS[4].SetStats(False)
139 self.hist_XYS[5] = ROOT.TH2F('XYS5', ' ;x;y', 10, -u2, -u1, 10, -u2, -u1)
140 self.hist_XYS[5].SetStats(False)
141 self.hist_XYS[6] = ROOT.TH2F('XYS6', ' ;x;y', 10, -u5, +u5, 10, -u4, -u3)
142 self.hist_XYS[6].SetStats(False)
143 self.hist_XYS[7] = ROOT.TH2F('XYS7', ' ;x;y', 10, +u1, +u2, 10, -u2, -u1)
144 self.hist_XYS[7].SetStats(False)
145
146 self.hist_ZYS = ROOT.TH2F('ZYS', ' ;z;y', 10, -150.0, 150.0, 10, 125.0, 425.0)
147 self.hist_ZYS.SetStats(False)
148 ROOT.gStyle.SetOptStat(10)
149
150 self.bklmXY = []
151 r0 = 201.9 + 0.5 * 4.4 # cm
152 dr = 9.1 # cm
153 tan0 = math.tan(math.pi / 8.0)
154 g = ROOT.TGraph()
155 g.SetPoint(0, -200.0, 0.0)
156 g.SetPoint(1, +200.0, 0.0)
157 g.SetLineColor(19)
158 g.SetLineWidth(1)
159 g.SetLineStyle(3)
160 self.bklmXY.append(g)
161 g = ROOT.TGraph()
162 g.SetPoint(0, 0.0, -200.0)
163 g.SetPoint(1, 0.0, +200.0)
164 g.SetLineColor(19)
165 g.SetLineWidth(1)
166 g.SetLineStyle(3)
167 self.bklmXY.append(g)
168 g = ROOT.TGraph()
169 g.SetPoint(0, -5.0, 0.0)
170 g.SetPoint(1, +5.0, 0.0)
171 g.SetPoint(2, 0.0, 0.0)
172 g.SetPoint(3, 0.0, +5.0)
173 g.SetPoint(4, 0.0, -5.0)
174 g.SetLineColor(1)
175 g.SetLineWidth(1)
176 self.bklmXY.append(g)
177 for layer in range(0, 15):
178 r = r0 + layer * dr
179 x = r * tan0
180 g = ROOT.TGraph()
181 g.SetPoint(0, +r, -x)
182 g.SetPoint(1, +r, +x)
183 g.SetPoint(2, +x, +r)
184 g.SetPoint(3, -x, +r)
185 g.SetPoint(4, -r, +x)
186 g.SetPoint(5, -r, -x)
187 g.SetPoint(6, -x, -r)
188 g.SetPoint(7, +x, -r)
189 g.SetPoint(8, +r, -x)
190 if layer < 2:
191 g.SetLineColor(18)
192 else:
193 g.SetLineColor(17)
194 if (layer % 5) == 0:
195 g.SetLineStyle(3)
196 g.SetLineWidth(1)
197 self.bklmXY.append(g)
198
199 self.bklmZY = []
200 rF = r0 + 14 * dr
201 x0 = r0 * tan0
202 z0 = 47.0 # cm
203 zL = 220.0 # cm
204 g = ROOT.TGraph()
205 g.SetPoint(0, -zL + z0 - 140.0, 0.0)
206 g.SetPoint(1, +zL + z0 + 70.0, 0.0)
207 g.SetLineColor(19)
208 g.SetLineWidth(1)
209 g.SetLineStyle(3)
210 self.bklmZY.append(g)
211 g = ROOT.TGraph()
212 g.SetPoint(0, 0.0, -315.0)
213 g.SetPoint(1, 0.0, +340.0)
214 g.SetLineColor(19)
215 g.SetLineWidth(1)
216 g.SetLineStyle(3)
217 self.bklmZY.append(g)
218 g = ROOT.TGraph()
219 g.SetPoint(0, -5.0, 0.0)
220 g.SetPoint(1, +5.0, 0.0)
221 g.SetPoint(2, 0.0, 0.0)
222 g.SetPoint(3, 0.0, +5.0)
223 g.SetPoint(4, 0.0, -5.0)
224 g.SetLineColor(1)
225 g.SetLineWidth(1)
226 self.bklmZY.append(g)
227 g = ROOT.TGraph()
228 g.SetPoint(0, -zL + z0, +x0)
229 g.SetPoint(1, -zL + z0, +r0)
230 g.SetLineColor(18)
231 g.SetLineWidth(1)
232 g.SetLineStyle(3)
233 self.bklmZY.append(g)
234 g = ROOT.TGraph()
235 g.SetPoint(0, -zL + z0, -x0)
236 g.SetPoint(1, -zL + z0, -r0)
237 g.SetLineColor(18)
238 g.SetLineWidth(1)
239 g.SetLineStyle(3)
240 self.bklmZY.append(g)
241 g = ROOT.TGraph()
242 g.SetPoint(0, +zL + z0, +x0)
243 g.SetPoint(1, +zL + z0, +r0)
244 g.SetLineColor(18)
245 g.SetLineWidth(1)
246 g.SetLineStyle(3)
247 self.bklmZY.append(g)
248 g = ROOT.TGraph()
249 g.SetPoint(0, +zL + z0, -x0)
250 g.SetPoint(1, +zL + z0, -r0)
251 g.SetLineColor(18)
252 g.SetLineWidth(1)
253 g.SetLineStyle(3)
254 self.bklmZY.append(g)
255 g = ROOT.TGraph()
256 g.SetPoint(0, -zL + z0, r0)
257 g.SetPoint(1, +zL + z0, r0)
258 g.SetPoint(2, +zL + z0, rF)
259 g.SetPoint(3, -zL + z0, rF)
260 g.SetPoint(4, -zL + z0, r0)
261 g.SetLineColor(18)
262 g.SetLineWidth(1)
263 self.bklmZY.append(g)
264 g = ROOT.TGraph()
265 g.SetPoint(0, -zL + z0, -r0)
266 g.SetPoint(1, +zL + z0, -r0)
267 g.SetPoint(2, +zL + z0, -rF)
268 g.SetPoint(3, -zL + z0, -rF)
269 g.SetPoint(4, -zL + z0, -r0)
270 g.SetLineColor(18)
271 g.SetLineWidth(1)
272 self.bklmZY.append(g)
273 g = ROOT.TGraph()
274 g.SetPoint(0, -zL + z0, -x0)
275 g.SetPoint(1, +zL + z0, -x0)
276 g.SetPoint(2, +zL + z0, +x0)
277 g.SetPoint(3, -zL + z0, +x0)
278 g.SetPoint(4, -zL + z0, -x0)
279 g.SetLineColor(18)
280 g.SetLineWidth(1)
281 self.bklmZY.append(g)
282
283 self.bklmZYL = []
284 for layer in range(0, 15):
285 r = r0 + layer * dr
286 z0 = 47.0 # cm
287 zL = 220.0 # cm
288 g = ROOT.TGraph()
289 g.SetPoint(0, -zL + z0, r)
290 g.SetPoint(1, +zL + z0, r)
291 g.SetLineColor(19)
292 g.SetLineWidth(1)
293 self.bklmZYL.append(g)
294
295 self.electIdToModuleId = bklmDB.fillDB()
296
297 self.sectorFBToDC = [11, 15, 2, 6, 10, 14, 3, 7, 9, 13, 0, 4, 8, 12, 1, 5]
298
299 self.dcToSectorFB = [10, 14, 2, 6, 11, 15, 3, 7, 12, 8, 4, 0, 13, 9, 5, 1]
300
301 self.t0Cal = 312
302
303 self.t0Cal2d = 308
304
305 self.ct0Cal = 455
306
307 self.ct0Cal2d = 520
308
309 self.t0RPC = [8, -14, -6, -14, -2, 10, 9, 13, 0, -10, -14, -20, 2, 6, 14, 11]
310
311 self.ct0Scint = [-1, -33, -46, -33, -2, 32, 51, 32, 0, -32, -45, -33, -4, 34, 45, 27]
312
def fillDB()
Definition: bklmDB.py:16

◆ terminate()

def terminate (   self)
Handle job termination: close event-display file

Definition at line 313 of file EventDisplayer.py.

313 def terminate(self):
314 """Handle job termination: close event-display file"""
315
316 pdfNameLast = f'{self.eventPdfName}]'
317 self.eventCanvas.Print(pdfNameLast, self.lastTitle)
318 print('Goodbye')
319

Member Data Documentation

◆ BKLM_ID

int BKLM_ID = 0x07000000
static

COPPER base identifier for BKLM readout.

Definition at line 32 of file EventDisplayer.py.

◆ BKLM_LAYER_BIT

int BKLM_LAYER_BIT = 7
static

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

Definition at line 40 of file EventDisplayer.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 52 of file EventDisplayer.py.

◆ BKLM_MAXSTRIP_BIT

int BKLM_MAXSTRIP_BIT = 15
static

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

Definition at line 46 of file EventDisplayer.py.

◆ BKLM_MAXSTRIP_MASK

tuple BKLM_MAXSTRIP_MASK = (63 << BKLM_MAXSTRIP_BIT)
static

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

Definition at line 58 of file EventDisplayer.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 60 of file EventDisplayer.py.

◆ BKLM_PLANE_BIT

int BKLM_PLANE_BIT = 6
static

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

Definition at line 38 of file EventDisplayer.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 50 of file EventDisplayer.py.

◆ BKLM_SECTION_BIT

int BKLM_SECTION_BIT = 14
static

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

Definition at line 44 of file EventDisplayer.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 56 of file EventDisplayer.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 42 of file EventDisplayer.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 54 of file EventDisplayer.py.

◆ BKLM_STRIP_BIT

int BKLM_STRIP_BIT = 0
static

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

Definition at line 36 of file EventDisplayer.py.

◆ BKLM_STRIP_MASK

int BKLM_STRIP_MASK = 0x3f
static

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

Definition at line 48 of file EventDisplayer.py.

◆ bklmXY

bklmXY

list of line-segment (x,y) points for the BKLM end view

Definition at line 150 of file EventDisplayer.py.

◆ bklmZY

bklmZY

list of line-segment (z,y) points for the BKLM side view

Definition at line 199 of file EventDisplayer.py.

◆ bklmZYL

bklmZYL

list of line-segment (z,y) points for the BKLM sector's zoomed and rotated side view

Definition at line 283 of file EventDisplayer.py.

◆ cosine

cosine

table of cosines for the BKLM-sector normals

Definition at line 108 of file EventDisplayer.py.

◆ ct0Cal

ct0Cal

scint-ctime calibration adjustment (ns) for rawKLMs

Definition at line 305 of file EventDisplayer.py.

◆ ct0Cal2d

ct0Cal2d

scint-ctime calibration adjustment (ns) for KLMHit2ds

Definition at line 307 of file EventDisplayer.py.

◆ ct0Scint

ct0Scint

per-sector variations in scint-ctime calibration adjustment (ns) for rawKLMs

Definition at line 311 of file EventDisplayer.py.

◆ dcToSectorFB

dcToSectorFB

map for data concentrator -> sectorFB

Definition at line 299 of file EventDisplayer.py.

◆ EKLM_ID

int EKLM_ID = 0x08000000
static

COPPER base identifier for EKLM readout.

Definition at line 34 of file EventDisplayer.py.

◆ electIdToModuleId

electIdToModuleId

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

Definition at line 295 of file EventDisplayer.py.

◆ eventCanvas

eventCanvas

TCanvas on which event displays will be drawn.

Definition at line 99 of file EventDisplayer.py.

◆ eventCounter

eventCounter

event counter (needed for PDF table of contents' ordinal event#)

Definition at line 87 of file EventDisplayer.py.

◆ eventDisplays

eventDisplays

event-display counter

Definition at line 89 of file EventDisplayer.py.

◆ eventPdfName

eventPdfName

internal copy of the pathname of the output event-display PDF file

Definition at line 79 of file EventDisplayer.py.

◆ exp

exp

internal copy of experiment number

Definition at line 75 of file EventDisplayer.py.

◆ hist_XY

hist_XY

blank scatterplot to define the bounds of the BKLM end view

Definition at line 116 of file EventDisplayer.py.

◆ hist_XYS

hist_XYS

list of blank scatterplots to define the per-sector bounds of the BKLM end view

Definition at line 128 of file EventDisplayer.py.

◆ hist_ZY

hist_ZY

blank scatterplot to define the bounds of the BKLM side view

Definition at line 119 of file EventDisplayer.py.

◆ hist_ZYS

hist_ZYS

blank scatterplot to define the per-sector bounds of the rotated BKLM side view

Definition at line 146 of file EventDisplayer.py.

◆ lastTitle

lastTitle

title of the last-drawn event display (needed for PDF table of contents' last event)

Definition at line 91 of file EventDisplayer.py.

◆ maxDisplays

maxDisplays

internal copy of the maximum number of event displays to write

Definition at line 81 of file EventDisplayer.py.

◆ minMuidHits

minMuidHits

internal copy of the minimum number of MuidHits in the event for event display

Definition at line 85 of file EventDisplayer.py.

◆ minRPCHits

minRPCHits

internal copy of the minimum number of RPC KLMHit2ds in any sector for event display

Definition at line 83 of file EventDisplayer.py.

◆ run

run

internal copy of run number

Definition at line 77 of file EventDisplayer.py.

◆ sectorFBToDC

sectorFBToDC

map for sectorFB -> data concentrator

Definition at line 297 of file EventDisplayer.py.

◆ sine

sine

table of sines for the BKLM-sector normals

Definition at line 110 of file EventDisplayer.py.

◆ t0Cal

t0Cal

RPC-time calibration adjustment (ns) for rawKLMs.

Definition at line 301 of file EventDisplayer.py.

◆ t0Cal2d

t0Cal2d

RPC-time calibration adjustment (ns) for KLMHit2ds.

Definition at line 303 of file EventDisplayer.py.

◆ t0RPC

t0RPC

per-sector variations in RPC-time calibration adjustment (ns) for rawKLMs

Definition at line 309 of file EventDisplayer.py.


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