Belle II Software  release-08-01-10
TOPTBCResolution Class Reference
Inheritance diagram for TOPTBCResolution:
Collaboration diagram for TOPTBCResolution:

Public Member Functions

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

Public Attributes

 outname
 output name
 
 m_calpulseMaxWidth
 output name
 
 m_calpulseMinWidth
 output name
 
 m_calpulseMaxAmp
 output name
 
 m_calpulseMinAmp
 output name
 
 m_ignoreNotCalibrated
 output name
 

Static Public Attributes

 h_WidthVSAmplitude_1
 Width VS amplitude, first calibration pulse. More...
 
 h_WidthVSAmplitude_2
 Width VS amplitude, second calibration pulse. More...
 
 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 More...
 
 h_dVdtRisingVSdVdtFalling_2
 dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse More...
 
 h_dVdtRisingDifference
 Difference between the dV/dt of the first and the second calpulse, using the rising edges. More...
 
 h_dVdtFallingDifference
 Difference between the dV/dt of the first and the second calpulse, using the rising edges. More...
 
 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. More...
 
 h_DeltaTVSChannel_FF
 DeltaT falling-falling VS channel. More...
 
 h_DeltaTVSChannel_FR
 DeltaT falling-rising VS channel. More...
 
 h_DeltaTVSChannel_RF
 DeltaT rising-falling VS channel. More...
 
 h_DeltaTVSdVdt_RR
 DeltaT rising-rising VS average of dV/dt on the first and second pulse. More...
 
 h_DeltaTVSdVdt_FF
 DeltaT falling-falling VS average of dV/dt on the first and second pulse. More...
 
 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)
 
string 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 27 of file studyTBCResolution.py.

Member Function Documentation

◆ event()

def event (   self)
 Event processor: fill histograms 

Definition at line 155 of file studyTBCResolution.py.

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

◆ ignoreNotCalibrated()

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

Definition at line 150 of file studyTBCResolution.py.

◆ setMaxAmp()

def setMaxAmp (   self,
  maxAmp 
)
 Sets the maximum calpulse amplitude 

Definition at line 140 of file studyTBCResolution.py.

◆ setMaxWidth()

def setMaxWidth (   self,
  maxWidth 
)
 Sets the maximum calpulse width 

Definition at line 130 of file studyTBCResolution.py.

◆ setMinAmp()

def setMinAmp (   self,
  minAmp 
)
 Sets the minimum calpulse amplitude 

Definition at line 145 of file studyTBCResolution.py.

◆ setMinWidth()

def setMinWidth (   self,
  minWidth 
)
 Sets the minimum calpulse width 

Definition at line 135 of file studyTBCResolution.py.

◆ setOutputName()

def setOutputName (   self,
  outputname 
)
 Sets the output file name 

Definition at line 125 of file studyTBCResolution.py.

◆ terminate()

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

Definition at line 227 of file studyTBCResolution.py.

Member Data Documentation

◆ 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 82 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 87 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 92 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 77 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 102 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 97 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 64 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 61 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 51 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 56 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 32 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 37 of file studyTBCResolution.py.


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