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

Public Member Functions

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

Static Public Attributes

 h_WidthVSAmplitude_1
 Width VS amplitude, first calibration pulse.
 
 h_WidthVSAmplitude_2
 Width VS amplitude, second calibration pulse.
 
 h_dVdtRising_1 = TH1F('dVdtRising_1', ' dV/dt of the TOPRawDigits (rising edge), first calibration pulse', 1000, 0, 1000)
 dV/dt on the rising edge, first calibration pulse
 
 h_dVdtRising_2 = TH1F('dVdtRising_2', ' dV/dt of the TOPRawDigits (rising edge), second calibration pulse', 1000, 0, 1000)
 dV/dt on the rising edge, second calibration pulse
 
 h_dVdtFalling_1 = TH1F('dVdtFalling_1', ' dV/dt of the TOPRawDigits (falling edge), first calibration pulse', 1000, 0, 1000)
 dV/dt on the falling edge, first calibration pulse
 
 h_dVdtFalling_2 = TH1F('dVdtFalling_2', ' dV/dt of the TOPRawDigits (falling edge), second calibration pulse', 1000, 0, 1000)
 dV/dt on the falling edge, second calibration pulse
 
 h_dVdtRisingVSdVdtFalling_1
 dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse
 
 h_dVdtRisingVSdVdtFalling_2
 dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse
 
 h_dVdtRisingDifference
 Difference between the dV/dt of the first and the second calpulse, using the rising edges.
 
 h_dVdtFallingDifference
 Difference between the dV/dt of the first and the second calpulse, using the rising edges.
 
 h_DeltaT_RR = TH1F('DeltaT_RR', ' DeltaT bewteen the rising edges', 4000, 10, 30)
 DeltaT rising-rising.
 
 h_DeltaT_FF = TH1F('DeltaT_FF', ' DeltaT bewteen the falling edges', 4000, 10, 30)
 DeltaT falling-falling.
 
 h_DeltaT_FR = TH1F('DeltaT_FR', ' DeltaT bewteen falling and rising edges', 4000, 10, 30)
 DeltaT falling-rising.
 
 h_DeltaT_RF = TH1F('DeltaT_RF', ' DeltaT bewteen rising and falling edges', 4000, 10, 30)
 DeltaT rising-falling.
 
 h_DeltaTVSChannel_RR
 DeltaT rising-rising VS channel.
 
 h_DeltaTVSChannel_FF
 DeltaT falling-falling VS channel.
 
 h_DeltaTVSChannel_FR
 DeltaT falling-rising VS channel.
 
 h_DeltaTVSChannel_RF
 DeltaT rising-falling VS channel.
 
 h_DeltaTVSdVdt_RR
 DeltaT rising-rising VS average of dV/dt on the first and second pulse.
 
 h_DeltaTVSdVdt_FF
 DeltaT falling-falling VS average of dV/dt on the first and second pulse.
 
 h_ResolutionVSdVdt_FF = TGraphErrors()
 DeltaT resolution VS average of dV/dt (falling-falling)
 
 h_ResolutionVSdVdt_RR = TGraphErrors()
 DeltaT resolution VS average of dV/dt (rising-rising)
 
str outname = 'outStudyTBCResolution.root'
 output root file
 
int m_calpulseMaxWidth = 3.
 maximum width to flag a calpulse candidate
 
float m_calpulseMinWidth = 0.5
 minimum width to flag a calpulse candidate
 
int m_calpulseMaxAmp = 700.
 minimum amplitude to flag a calpulse candidate
 
int m_calpulseMinAmp = 250.
 minimum amplitude to flag a calpulse candidate
 
bool m_ignoreNotCalibrated = True
 ignores the hits wthout calibration
 

Detailed Description

 Module to study resolution and performances of the TOP Time Base Calibration.

Definition at line 26 of file studyTBCResolution.py.

Member Function Documentation

◆ event()

event ( self)
 Event processor: fill histograms 

Definition at line 154 of file studyTBCResolution.py.

154 def event(self):
155 ''' Event processor: fill histograms '''
156
157 digits = Belle2.PyStoreArray('TOPDigits')
158
159 for ipulse1, digit1 in enumerate(digits):
160 if (self.ignoreNotCalibrated and not digit1.isTimeBaseCalibrated()):
161 continue
162 if (digit1.getHitQuality() != 0 and
163 digit1.getPulseHeight() > self.m_calpulseMinAmp and
164 digit1.getPulseHeight() < self.m_calpulseMaxAmp and
165 digit1.getPulseWidth() > self.m_calpulseMinWidth and
166 digit1.getPulseWidth() < self.m_calpulseMaxWidth):
167
168 slotID = digit1.getModuleID()
169 hwchan = digit1.getChannel()
170 for ipulse2, digit2 in enumerate(digits, start=ipulse1 + 1):
171 if (self.ignoreNotCalibrated and not digit2.isTimeBaseCalibrated()):
172 continue
173
174 if (digit2.getHitQuality() != 0 and
175 digit2.getPulseHeight() > self.m_calpulseMinAmp and
176 digit2.getPulseHeight() < self.m_calpulseMaxAmp and
177 digit2.getPulseWidth() > self.m_calpulseMinWidth and
178 digit2.getPulseWidth() < self.m_calpulseMaxWidth and
179 slotID == digit2.getModuleID() and
180 hwchan == digit2.getChannel()):
181
182 # finds which one is the first calpulse
183 rawDigitFirst = digit1.getRelated('TOPRawDigits')
184 rawDigitSecond = digit2.getRelated('TOPRawDigits')
185 if digit1.getTime() > digit2.getTime():
186 rawDigitFirst = digit2.getRelated('TOPRawDigits')
187 rawDigitSecond = digit1.getRelated('TOPRawDigits')
188 digitFirst = rawDigitFirst.getRelated('TOPDigits')
189 digitSecond = rawDigitSecond.getRelated('TOPDigits')
190
191 globalCh = hwchan + 512 * (slotID - 1)
192 dV1_R = rawDigitFirst.getLeadingSlope()
193 dV1_F = -rawDigitFirst.getFallingSlope()
194 dV2_R = rawDigitSecond.getLeadingSlope()
195 dV2_F = -rawDigitSecond.getFallingSlope()
196 t1_R = digitFirst.getTime()
197 t1_F = digitFirst.getTime() + digitFirst.getPulseWidth()
198 t2_R = digitSecond.getTime()
199 t2_F = digitSecond.getTime() + digitSecond.getPulseWidth()
200 amp1 = digitFirst.getPulseHeight()
201 amp2 = digitSecond.getPulseHeight()
202 w1 = digitFirst.getPulseWidth()
203 w2 = digitSecond.getPulseWidth()
204
205 self.h_WidthVSAmplitude_1.Fill(amp1, w1)
206 self.h_WidthVSAmplitude_2.Fill(amp2, w2)
207 self.h_dVdtRising_1.Fill(dV1_R)
208 self.h_dVdtRising_2.Fill(dV2_R)
209 self.h_dVdtFalling_1.Fill(dV1_F)
210 self.h_dVdtFalling_2.Fill(dV2_F)
211 self.h_dVdtRisingVSdVdtFalling_1.Fill(dV1_F, dV1_R)
212 self.h_dVdtRisingVSdVdtFalling_2.Fill(dV2_F, dV2_R)
213 self.h_dVdtRisingDifference.Fill(dV1_R - dV2_R)
214 self.h_dVdtFallingDifference.Fill(dV1_F - dV2_F)
215 self.h_DeltaT_RR.Fill(t2_R - t1_R)
216 self.h_DeltaT_FF.Fill(t2_F - t1_F)
217 self.h_DeltaT_FR.Fill(t2_F - t1_R)
218 self.h_DeltaT_RF.Fill(t2_R - t1_F)
219 self.h_DeltaTVSChannel_RR.Fill(globalCh, t2_R - t1_R)
220 self.h_DeltaTVSChannel_FF.Fill(globalCh, t2_F - t1_F)
221 self.h_DeltaTVSChannel_FR.Fill(globalCh, t2_F - t1_R)
222 self.h_DeltaTVSChannel_RF.Fill(globalCh, t2_R - t1_F)
223 self.h_DeltaTVSdVdt_RR.Fill(0.5 * (dV1_R + dV2_R), t2_R - t1_R)
224 self.h_DeltaTVSdVdt_FF.Fill(0.5 * (dV1_F + dV2_F), t2_F - t1_F)
225
A (simplified) python wrapper for StoreArray.

◆ ignoreNotCalibrated()

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

Definition at line 149 of file studyTBCResolution.py.

149 def ignoreNotCalibrated(self, ignoreNotCal):
150 ''' Sets the flag to ingore the hits without calibration '''
151
152 self.m_ignoreNotCalibrated = ignoreNotCal
153

◆ setMaxAmp()

setMaxAmp ( self,
maxAmp )
 Sets the maximum calpulse amplitude 

Definition at line 139 of file studyTBCResolution.py.

139 def setMaxAmp(self, maxAmp):
140 ''' Sets the maximum calpulse amplitude '''
141
142 self.m_calpulseMaxAmp = maxAmp
143

◆ setMaxWidth()

setMaxWidth ( self,
maxWidth )
 Sets the maximum calpulse width 

Definition at line 129 of file studyTBCResolution.py.

129 def setMaxWidth(self, maxWidth):
130 ''' Sets the maximum calpulse width '''
131
132 self.m_calpulseMaxWidth = maxWidth
133

◆ setMinAmp()

setMinAmp ( self,
minAmp )
 Sets the minimum calpulse amplitude 

Definition at line 144 of file studyTBCResolution.py.

144 def setMinAmp(self, minAmp):
145 ''' Sets the minimum calpulse amplitude '''
146
147 self.m_calpulseMinAmp = minAmp
148

◆ setMinWidth()

setMinWidth ( self,
minWidth )
 Sets the minimum calpulse width 

Definition at line 134 of file studyTBCResolution.py.

134 def setMinWidth(self, minWidth):
135 ''' Sets the minimum calpulse width '''
136
137 self.m_calpulseMinWidth = minWidth
138

◆ setOutputName()

setOutputName ( self,
outputname )
 Sets the output file name 

Definition at line 124 of file studyTBCResolution.py.

124 def setOutputName(self, outputname):
125 ''' Sets the output file name '''
126
127 self.outname = outputname
128

◆ terminate()

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

Definition at line 226 of file studyTBCResolution.py.

226 def terminate(self):
227 ''' Write histograms to file, fills and fits the resolution plots'''
228
229 self.h_ResolutionVSdVdt_RR.SetName('ResolutionVSdVdt_RR')
230 self.h_ResolutionVSdVdt_RR.SetTitle('Resolution VS dV/dt (rising-rising)')
231 self.h_ResolutionVSdVdt_FF.SetName('ResolutionVSdVdt_FF')
232 self.h_ResolutionVSdVdt_FF.SetTitle('Resolution VS dV/dt (falling-falling)')
233
234 for ibin in range(0, 10):
235 projection = self.h_DeltaTVSdVdt_RR.ProjectionY("tmpProj", ibin * 100 + 1, (ibin + 1) * 100)
236 gaussFit = TF1("gaussFit", "[0]*exp(-0.5*((x-[1])/[2])**2)", 10., 30.)
237 gaussFit.SetParameter(0, 1.)
238 gaussFit.SetParameter(1, projection.GetMean())
239 gaussFit.SetParameter(2, projection.GetRMS())
240 gaussFit.SetParLimits(2, 0., 3. * projection.GetRMS())
241
242 projection.Fit("gaussFit")
243 self.h_ResolutionVSdVdt_RR.SetPoint(ibin, ibin * 100. + 50., gaussFit.GetParameter(2))
244 self.h_ResolutionVSdVdt_RR.SetPointError(ibin, 50., gaussFit.GetParError(2))
245
246 tfile = TFile(self.outname, 'recreate')
247 self.h_WidthVSAmplitude_1.GetXaxis().SetTitle("TOPDigit amplitude [ADC counts]")
248 self.h_WidthVSAmplitude_1.GetYaxis().SetTitle("TOPDigit width [ns]")
249 self.h_WidthVSAmplitude_1.Write()
250 self.h_WidthVSAmplitude_2.GetXaxis().SetTitle("TOPDigit amplitude [ADC counts]")
251 self.h_WidthVSAmplitude_2.GetYaxis().SetTitle("TOPDigit width [ns]")
252 self.h_WidthVSAmplitude_2.Write()
253
254 self.h_dVdtRising_1.GetXaxis().SetTitle("dV/dt [ADC counts / sample]")
255 self.h_dVdtRising_1.Write()
256 self.h_dVdtRising_2.GetXaxis().SetTitle("dV/dt [ADC counts / sample]")
257 self.h_dVdtRising_2.Write()
258 self.h_dVdtFalling_1.GetXaxis().SetTitle("dV/dt [ADC counts / sample]")
259 self.h_dVdtFalling_1.Write()
260 self.h_dVdtFalling_2.GetXaxis().SetTitle("dV/dt [ADC counts / sample]")
261 self.h_dVdtFalling_2.Write()
262
263 self.h_dVdtRisingVSdVdtFalling_1.GetXaxis().SetTitle("dV/dt on the falling edge [ADC counts / sample]")
264 self.h_dVdtRisingVSdVdtFalling_1.GetYaxis().SetTitle("dV/dt on the rising edge [ADC counts / sample]")
265 self.h_dVdtRisingVSdVdtFalling_1.Write()
266
267 self.h_dVdtRisingVSdVdtFalling_2.GetXaxis().SetTitle("dV/dt on the falling edge [ADC counts / sample]")
268 self.h_dVdtRisingVSdVdtFalling_2.GetYaxis().SetTitle("dV/dt on the rising edge [ADC counts / sample]")
269 self.h_dVdtRisingVSdVdtFalling_2.Write()
270
271 self.h_dVdtFallingDifference.GetXaxis().SetTitle("dV/dt_1 - dV/dt_2 [ADC counts / sample]")
272 self.h_dVdtFallingDifference.Write()
273
274 self.h_dVdtRisingDifference.GetXaxis().SetTitle("dV/dt_1 - dV/dt_2 [ADC counts / sample]")
275 self.h_dVdtRisingDifference.Write()
276
277 self.h_DeltaT_RR.GetXaxis().SetTitle("#Delta t [ns]")
278 self.h_DeltaT_RR.Write()
279 self.h_DeltaT_FF.GetXaxis().SetTitle("#Delta t [ns]")
280 self.h_DeltaT_FF.Write()
281 self.h_DeltaT_FR.GetXaxis().SetTitle("#Delta t [ns]")
282 self.h_DeltaT_FR.Write()
283 self.h_DeltaT_RF.GetXaxis().SetTitle("#Delta t [ns]")
284 self.h_DeltaT_RF.Write()
285
286 self.h_DeltaTVSChannel_RR.GetXaxis().SetTitle("Global channel number [hwChannel + 512*(slotID-1)]")
287 self.h_DeltaTVSChannel_RR.GetYaxis().SetTitle("#Delta t [ns]")
288 self.h_DeltaTVSChannel_RR.Write()
289 self.h_DeltaTVSChannel_FF.GetXaxis().SetTitle("Global channel number [hwChannel + 512*(slotID-1)]")
290 self.h_DeltaTVSChannel_FF.GetYaxis().SetTitle("#Delta t [ns]")
291 self.h_DeltaTVSChannel_FF.Write()
292 self.h_DeltaTVSChannel_FR.GetXaxis().SetTitle("Global channel number [hwChannel + 512*(slotID-1)]")
293 self.h_DeltaTVSChannel_FR.GetYaxis().SetTitle("#Delta t [ns]")
294 self.h_DeltaTVSChannel_FR.Write()
295 self.h_DeltaTVSChannel_RF.GetXaxis().SetTitle("Global channel number [hwChannel + 512*(slotID-1)]")
296 self.h_DeltaTVSChannel_RF.GetYaxis().SetTitle("#Delta t [ns]")
297 self.h_DeltaTVSChannel_RF.Write()
298
299 self.h_DeltaTVSdVdt_RR.GetXaxis().SetTitle("Average of dV/dt on first and second pulse [ACD counts / sample]")
300 self.h_DeltaTVSdVdt_RR.GetYaxis().SetTitle("#Delta t [ns]")
301 self.h_DeltaTVSdVdt_RR.Write()
302
303 self.h_DeltaTVSdVdt_FF.GetXaxis().SetTitle("Average of dV/dt on first and second pulse [ACD counts / sample]")
304 self.h_DeltaTVSdVdt_FF.GetYaxis().SetTitle("#Delta t [ns]")
305 self.h_DeltaTVSdVdt_FF.Write()
306
307 self.h_ResolutionVSdVdt_RR.Write()
308 tfile.Close()
309
310

Member Data Documentation

◆ h_DeltaT_FF

h_DeltaT_FF = TH1F('DeltaT_FF', ' DeltaT bewteen the falling edges', 4000, 10, 30)
static

DeltaT falling-falling.

Definition at line 69 of file studyTBCResolution.py.

◆ h_DeltaT_FR

h_DeltaT_FR = TH1F('DeltaT_FR', ' DeltaT bewteen falling and rising edges', 4000, 10, 30)
static

DeltaT falling-rising.

Definition at line 71 of file studyTBCResolution.py.

◆ h_DeltaT_RF

h_DeltaT_RF = TH1F('DeltaT_RF', ' DeltaT bewteen rising and falling edges', 4000, 10, 30)
static

DeltaT rising-falling.

Definition at line 73 of file studyTBCResolution.py.

◆ h_DeltaT_RR

h_DeltaT_RR = TH1F('DeltaT_RR', ' DeltaT bewteen the rising edges', 4000, 10, 30)
static

DeltaT rising-rising.

Definition at line 67 of file studyTBCResolution.py.

◆ h_DeltaTVSChannel_FF

h_DeltaTVSChannel_FF
static
Initial value:
= TH2F(
'DeltaTVSChannel_FF',
' DeltaT bewteen the falling edges, as function of the channel number',
512 * 16, 0, 512 * 16, 4000, 10., 30.)

DeltaT falling-falling VS channel.

Definition at line 81 of file studyTBCResolution.py.

◆ h_DeltaTVSChannel_FR

h_DeltaTVSChannel_FR
static
Initial value:
= TH2F(
'DeltaTVSChannel_FR',
' DeltaT bewteen falling (pulse 1) and rising (pulse 2) edge, as function of the channel number',
512 * 16, 0, 512 * 16, 4000, 10., 30.)

DeltaT falling-rising VS channel.

Definition at line 86 of file studyTBCResolution.py.

◆ h_DeltaTVSChannel_RF

h_DeltaTVSChannel_RF
static
Initial value:
= TH2F(
'DeltaTVSChannel_RF',
' DeltaT bewteen rising (pulse 1) and falling (pulse 2) edge, as function of the channel number',
512 * 16, 0, 512 * 16, 4000, 10., 30.)

DeltaT rising-falling VS channel.

Definition at line 91 of file studyTBCResolution.py.

◆ h_DeltaTVSChannel_RR

h_DeltaTVSChannel_RR
static
Initial value:
= TH2F(
'DeltaTVSChannel_RR',
' DeltaT bewteen the rising edges, as function of the channel number',
512 * 16, 0, 512 * 16, 4000, 10., 30.)

DeltaT rising-rising VS channel.

Definition at line 76 of file studyTBCResolution.py.

◆ h_DeltaTVSdVdt_FF

h_DeltaTVSdVdt_FF
static
Initial value:
= TH2F(
'DeltaTVSdVdt_FF',
'DeltaT bewteen the rising edges VS average of dV/dt on the first and second pulser',
1000, 0., 1000., 4000, 10., 30.)

DeltaT falling-falling VS average of dV/dt on the first and second pulse.

Definition at line 101 of file studyTBCResolution.py.

◆ h_DeltaTVSdVdt_RR

h_DeltaTVSdVdt_RR
static
Initial value:
= TH2F(
'DeltaTVSdVdt_RR',
'DeltaT bewteen the rising edges VS average of dV/dt on the first and second pulser',
1000, 0., 1000., 4000, 10., 30.)

DeltaT rising-rising VS average of dV/dt on the first and second pulse.

Definition at line 96 of file studyTBCResolution.py.

◆ h_dVdtFalling_1

h_dVdtFalling_1 = TH1F('dVdtFalling_1', ' dV/dt of the TOPRawDigits (falling edge), first calibration pulse', 1000, 0, 1000)
static

dV/dt on the falling edge, first calibration pulse

Definition at line 46 of file studyTBCResolution.py.

◆ h_dVdtFalling_2

h_dVdtFalling_2 = TH1F('dVdtFalling_2', ' dV/dt of the TOPRawDigits (falling edge), second calibration pulse', 1000, 0, 1000)
static

dV/dt on the falling edge, second calibration pulse

Definition at line 48 of file studyTBCResolution.py.

◆ h_dVdtFallingDifference

h_dVdtFallingDifference
static
Initial value:
= TH1F(
'dVdtFallingDifference', ' difference between the falling edge dV/dt of the first and the second pulse', 1000, -500, 500)

Difference between the dV/dt of the first and the second calpulse, using the rising edges.

Definition at line 63 of file studyTBCResolution.py.

◆ h_dVdtRising_1

h_dVdtRising_1 = TH1F('dVdtRising_1', ' dV/dt of the TOPRawDigits (rising edge), first calibration pulse', 1000, 0, 1000)
static

dV/dt on the rising edge, first calibration pulse

Definition at line 42 of file studyTBCResolution.py.

◆ h_dVdtRising_2

h_dVdtRising_2 = TH1F('dVdtRising_2', ' dV/dt of the TOPRawDigits (rising edge), second calibration pulse', 1000, 0, 1000)
static

dV/dt on the rising edge, second calibration pulse

Definition at line 44 of file studyTBCResolution.py.

◆ h_dVdtRisingDifference

h_dVdtRisingDifference
static
Initial value:
= TH1F(
'dVdtRisingDifference', ' difference between the rising edge dV/dt of the first and the second pulse', 1000, -500, 500)

Difference between the dV/dt of the first and the second calpulse, using the rising edges.

Definition at line 60 of file studyTBCResolution.py.

◆ h_dVdtRisingVSdVdtFalling_1

h_dVdtRisingVSdVdtFalling_1
static
Initial value:
= TH2F(
'dVdtRisingVSdVdtFalling_1',
' dV/dt of the TOPRawDigit: rising edge VS falling edge, first calibration pulse ',
1000, 0, 1000, 1000, 0., 1000.)

dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse

Definition at line 50 of file studyTBCResolution.py.

◆ h_dVdtRisingVSdVdtFalling_2

h_dVdtRisingVSdVdtFalling_2
static
Initial value:
= TH2F(
'dVdtRisingVSdVdtFalling_2',
' dV/dt of the TOPRawDigit: rising edge VS falling edge, second calibration pulse ',
1000, 0, 1000, 1000, 0., 1000.)

dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse

Definition at line 55 of file studyTBCResolution.py.

◆ h_ResolutionVSdVdt_FF

h_ResolutionVSdVdt_FF = TGraphErrors()
static

DeltaT resolution VS average of dV/dt (falling-falling)

Definition at line 107 of file studyTBCResolution.py.

◆ h_ResolutionVSdVdt_RR

h_ResolutionVSdVdt_RR = TGraphErrors()
static

DeltaT resolution VS average of dV/dt (rising-rising)

Definition at line 109 of file studyTBCResolution.py.

◆ h_WidthVSAmplitude_1

h_WidthVSAmplitude_1
static
Initial value:
= TH2F(
'WidthVSAmplitude_1',
'Width VS amplidute of the TOPDigits, first calibration pulse',
2000, 0., 2000, 100, 0., 10.)

Width VS amplitude, first calibration pulse.

Definition at line 31 of file studyTBCResolution.py.

◆ h_WidthVSAmplitude_2

h_WidthVSAmplitude_2
static
Initial value:
= TH2F(
'WidthVSAmplitude_2',
'Width VS amplidute of the TOPDigits, second calibration pulse',
2000, 0., 2000, 100, 0., 10.)

Width VS amplitude, second calibration pulse.

Definition at line 36 of file studyTBCResolution.py.

◆ m_calpulseMaxAmp

int m_calpulseMaxAmp = 700.
static

minimum amplitude to flag a calpulse candidate

Definition at line 118 of file studyTBCResolution.py.

◆ m_calpulseMaxWidth

m_calpulseMaxWidth = 3.
static

maximum width to flag a calpulse candidate

Definition at line 114 of file studyTBCResolution.py.

◆ m_calpulseMinAmp

int m_calpulseMinAmp = 250.
static

minimum amplitude to flag a calpulse candidate

Definition at line 120 of file studyTBCResolution.py.

◆ m_calpulseMinWidth

float m_calpulseMinWidth = 0.5
static

minimum width to flag a calpulse candidate

Definition at line 116 of file studyTBCResolution.py.

◆ m_ignoreNotCalibrated

bool m_ignoreNotCalibrated = True
static

ignores the hits wthout calibration

Definition at line 122 of file studyTBCResolution.py.

◆ outname

str outname = 'outStudyTBCResolution.root'
static

output root file

Definition at line 112 of file studyTBCResolution.py.


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