Belle II Software prerelease-10-00-00a
TOPLaserHistogrammerModule Class Reference
Inheritance diagram for TOPLaserHistogrammerModule:
Collaboration diagram for TOPLaserHistogrammerModule:

Public Member Functions

 setOutputName (self, outputname)
 
 setMaxWidth (self, maxWidth)
 
 setMinWidth (self, minWidth)
 
 setMaxAmp (self, maxAmp)
 
 setMinAmp (self, minAmp)
 
 setMCCorrectionsFile (self, MCfile)
 
 ignoreNotCalibrated (self, ignoreNotCal)
 
 event (self)
 
 terminate (self)
 

Static Public Attributes

 h_LaserTimingVSChannel
 Width as function of the sample number in each channel.
 
 h_LaserTimingVSChannelOneSlot
 Laser timing in as function of the channel number.
 
list h_crossOccupancy
 cross occupancy
 
str outname = 'outStudyLaserResolution.root'
 output root file
 
bool m_ignoreNotCalibrated = True
 ignores the hits wthout calibration
 
bool m_runOnData = True
 ignores the hits wthout calibration
 
int m_maxWidth = 3.
 maximum width to accept a TOPDigit
 
float m_minWidth = 0.5
 minimum width to accept a TOPDigit
 
int m_maxAmp = 700.
 minimum amplitude to accept a TOPDigit
 
int m_minAmp = 250.
 minimum amplitude to accept a TOPDigit
 
str m_mcCorrectionsFile = '/group/belle2/group/detector/TOP/calibration/MCreferences/t0MC.root'
 root file with the MC distribution of the laser light, to get the light path corrections
 
list m_MCPeaks = [[]]
 positions of the first and second peak
 

Detailed Description

 Module to study resolution and performance of the top laser calibration.

Definition at line 29 of file studyLaserLight.py.

Member Function Documentation

◆ event()

event ( self)
 Event processor: fill histograms 

Definition at line 121 of file studyLaserLight.py.

121 def event(self):
122 ''' Event processor: fill histograms '''
123
124 digits = Belle2.PyStoreArray('TOPDigits')
125 nhits = [0 for i in range(16)]
126 for digit in digits:
127 if(not self.ignoreNotCalibrated and not digit.isTimeBaseCalibrated()):
128 continue
129 if (digit.getHitQuality() == 1 and
130 digit.getPulseWidth() > self.m_minWidth and digit.getPulseWidth() < self.m_maxWidth and
131 digit.getPulseHeight() > self.m_minAmp and digit.getPulseHeight() < self.m_maxAmp):
132 slotID = digit.getModuleID()
133 hwchan = digit.getChannel()
134 self.h_LaserTimingVSChannel.Fill(512 * (slotID - 1) + hwchan, digit.getTime())
135 simhits = digit.getRelationsWith('TOPSimHits')
136 nhits[slotID - 1] = nhits[slotID - 1] + 1
137 for simhit in simhits:
138 self.h_LaserTimingVSChannelOneSlot.Fill(hwchan, simhit.getTime())
139 for slotA in range(16):
140 for slotB in range(16):
141 self.h_crossOccupancy[slotA][slotB].Fill(nhits[slotA], nhits[slotB])
142
A (simplified) python wrapper for StoreArray.

◆ ignoreNotCalibrated()

ignoreNotCalibrated ( self,
ignoreNotCal )
 Sets the flag to ingore the hits without calibration 

Definition at line 116 of file studyLaserLight.py.

116 def ignoreNotCalibrated(self, ignoreNotCal):
117 ''' Sets the flag to ingore the hits without calibration '''
118
119 self.m_ignoreNotCalibrated = ignoreNotCal
120

◆ setMaxAmp()

setMaxAmp ( self,
maxAmp )
 Sets the maximum calpulse amplitude 

Definition at line 101 of file studyLaserLight.py.

101 def setMaxAmp(self, maxAmp):
102 ''' Sets the maximum calpulse amplitude '''
103
104 self.m_maxAmp = maxAmp
105

◆ setMaxWidth()

setMaxWidth ( self,
maxWidth )
 Sets the maximum calpulse width 

Definition at line 91 of file studyLaserLight.py.

91 def setMaxWidth(self, maxWidth):
92 ''' Sets the maximum calpulse width '''
93
94 self.m_maxWidth = maxWidth
95

◆ setMCCorrectionsFile()

setMCCorrectionsFile ( self,
MCfile )
 Sets the file containing the MC correction

Definition at line 111 of file studyLaserLight.py.

111 def setMCCorrectionsFile(self, MCfile):
112 ''' Sets the file containing the MC correction'''
113
114 self.m_mcCorrectionsFile = MCfile
115

◆ setMinAmp()

setMinAmp ( self,
minAmp )
 Sets the minimum calpulse amplitude 

Definition at line 106 of file studyLaserLight.py.

106 def setMinAmp(self, minAmp):
107 ''' Sets the minimum calpulse amplitude '''
108
109 self.m_minAmp = minAmp
110

◆ setMinWidth()

setMinWidth ( self,
minWidth )
 Sets the minimum calpulse width 

Definition at line 96 of file studyLaserLight.py.

96 def setMinWidth(self, minWidth):
97 ''' Sets the minimum calpulse width '''
98
99 self.m_minWidth = minWidth
100

◆ setOutputName()

setOutputName ( self,
outputname )
 Sets the output file name 

Definition at line 86 of file studyLaserLight.py.

86 def setOutputName(self, outputname):
87 ''' Sets the output file name '''
88
89 self.outname = outputname
90

◆ terminate()

terminate ( self)
 Write histograms to file, fills and fits the resolution plots

Definition at line 143 of file studyLaserLight.py.

143 def terminate(self):
144 ''' Write histograms to file, fills and fits the resolution plots'''
145 tfile = TFile(self.outname, 'recreate')
146 self.h_LaserTimingVSChannel.Write()
147 self.h_LaserTimingVSChannelOneSlot.Write()
148 for slotA in range(16):
149 for slotB in range(16):
150 self.h_crossOccupancy[slotA][slotB].Write()
151 tfile.Close()
152
153

Member Data Documentation

◆ h_crossOccupancy

list h_crossOccupancy
static
Initial value:
= [[TH2F(
'crossOccupancy_' + str(slotA) + '_' + str(slotB),
' ',
200,
0,
200,
200,
0.,
200) for slotA in range(16)] for slotB in range(16)]

cross occupancy

Definition at line 55 of file studyLaserLight.py.

◆ h_LaserTimingVSChannel

h_LaserTimingVSChannel
static
Initial value:
= TH2F(
'LaserTimingVSChannel',
'Laser timing in as function of the channel number',
512 * 16,
0,
512 * 16,
10000,
0.,
100)

Width as function of the sample number in each channel.

Definition at line 34 of file studyLaserLight.py.

◆ h_LaserTimingVSChannelOneSlot

h_LaserTimingVSChannelOneSlot
static
Initial value:
= TH2F(
'LaserTimingVSChannelOneSlot',
'Laser timing in as function of the channel number',
512,
0,
512,
10000,
0.,
100)

Laser timing in as function of the channel number.

Definition at line 44 of file studyLaserLight.py.

◆ m_ignoreNotCalibrated

bool m_ignoreNotCalibrated = True
static

ignores the hits wthout calibration

Definition at line 69 of file studyLaserLight.py.

◆ m_maxAmp

m_maxAmp = 700.
static

minimum amplitude to accept a TOPDigit

Definition at line 78 of file studyLaserLight.py.

◆ m_maxWidth

int m_maxWidth = 3.
static

maximum width to accept a TOPDigit

Definition at line 74 of file studyLaserLight.py.

◆ m_mcCorrectionsFile

str m_mcCorrectionsFile = '/group/belle2/group/detector/TOP/calibration/MCreferences/t0MC.root'
static

root file with the MC distribution of the laser light, to get the light path corrections

Definition at line 82 of file studyLaserLight.py.

◆ m_MCPeaks

list m_MCPeaks = [[]]
static

positions of the first and second peak

Definition at line 84 of file studyLaserLight.py.

◆ m_minAmp

int m_minAmp = 250.
static

minimum amplitude to accept a TOPDigit

Definition at line 80 of file studyLaserLight.py.

◆ m_minWidth

float m_minWidth = 0.5
static

minimum width to accept a TOPDigit

Definition at line 76 of file studyLaserLight.py.

◆ m_runOnData

bool m_runOnData = True
static

ignores the hits wthout calibration

Definition at line 71 of file studyLaserLight.py.

◆ outname

str outname = 'outStudyLaserResolution.root'
static

output root file

Definition at line 66 of file studyLaserLight.py.


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