22from ROOT
import Belle2
23from ROOT
import TH1F, TH2F, TF1, TFile, TGraphErrors
28 ''' Module to study resolution and performances of the TOP Time Base Calibration.'''
31 h_WidthVSAmplitude_1 = TH2F(
33 'Width VS amplidute of the TOPDigits, first calibration pulse',
34 2000, 0., 2000, 100, 0., 10.)
36 h_WidthVSAmplitude_2 = TH2F(
38 'Width VS amplidute of the TOPDigits, second calibration pulse',
39 2000, 0., 2000, 100, 0., 10.)
42 h_dVdtRising_1 = TH1F(
'dVdtRising_1',
' dV/dt of the TOPRawDigits (rising edge), first calibration pulse', 1000, 0, 1000)
44 h_dVdtRising_2 = TH1F(
'dVdtRising_2',
' dV/dt of the TOPRawDigits (rising edge), second calibration pulse', 1000, 0, 1000)
46 h_dVdtFalling_1 = TH1F(
'dVdtFalling_1',
' dV/dt of the TOPRawDigits (falling edge), first calibration pulse', 1000, 0, 1000)
48 h_dVdtFalling_2 = TH1F(
'dVdtFalling_2',
' dV/dt of the TOPRawDigits (falling edge), second calibration pulse', 1000, 0, 1000)
50 h_dVdtRisingVSdVdtFalling_1 = TH2F(
51 'dVdtRisingVSdVdtFalling_1',
52 ' dV/dt of the TOPRawDigit: rising edge VS falling edge, first calibration pulse ',
53 1000, 0, 1000, 1000, 0., 1000.)
55 h_dVdtRisingVSdVdtFalling_2 = TH2F(
56 'dVdtRisingVSdVdtFalling_2',
57 ' dV/dt of the TOPRawDigit: rising edge VS falling edge, second calibration pulse ',
58 1000, 0, 1000, 1000, 0., 1000.)
60 h_dVdtRisingDifference = TH1F(
61 'dVdtRisingDifference',
' difference between the rising edge dV/dt of the first and the second pulse', 1000, -500, 500)
63 h_dVdtFallingDifference = TH1F(
64 'dVdtFallingDifference',
' difference between the falling edge dV/dt of the first and the second pulse', 1000, -500, 500)
67 h_DeltaT_RR = TH1F(
'DeltaT_RR',
' DeltaT bewteen the rising edges', 4000, 10, 30)
69 h_DeltaT_FF = TH1F(
'DeltaT_FF',
' DeltaT bewteen the falling edges', 4000, 10, 30)
71 h_DeltaT_FR = TH1F(
'DeltaT_FR',
' DeltaT bewteen falling and rising edges', 4000, 10, 30)
73 h_DeltaT_RF = TH1F(
'DeltaT_RF',
' DeltaT bewteen rising and falling edges', 4000, 10, 30)
76 h_DeltaTVSChannel_RR = TH2F(
78 ' DeltaT bewteen the rising edges, as function of the channel number',
79 512 * 16, 0, 512 * 16, 4000, 10., 30.)
81 h_DeltaTVSChannel_FF = TH2F(
83 ' DeltaT bewteen the falling edges, as function of the channel number',
84 512 * 16, 0, 512 * 16, 4000, 10., 30.)
86 h_DeltaTVSChannel_FR = TH2F(
88 ' DeltaT bewteen falling (pulse 1) and rising (pulse 2) edge, as function of the channel number',
89 512 * 16, 0, 512 * 16, 4000, 10., 30.)
91 h_DeltaTVSChannel_RF = TH2F(
93 ' DeltaT bewteen rising (pulse 1) and falling (pulse 2) edge, as function of the channel number',
94 512 * 16, 0, 512 * 16, 4000, 10., 30.)
96 h_DeltaTVSdVdt_RR = TH2F(
98 'DeltaT bewteen the rising edges VS average of dV/dt on the first and second pulser',
99 1000, 0., 1000., 4000, 10., 30.)
101 h_DeltaTVSdVdt_FF = TH2F(
103 'DeltaT bewteen the rising edges VS average of dV/dt on the first and second pulser',
104 1000, 0., 1000., 4000, 10., 30.)
107 h_ResolutionVSdVdt_FF = TGraphErrors()
109 h_ResolutionVSdVdt_RR = TGraphErrors()
112 outname =
'outStudyTBCResolution.root'
114 m_calpulseMaxWidth = 3.
116 m_calpulseMinWidth = 0.5
118 m_calpulseMaxAmp = 700.
120 m_calpulseMinAmp = 250.
122 m_ignoreNotCalibrated =
True
125 ''' Sets the output file name '''
130 ''' Sets the maximum calpulse width '''
135 ''' Sets the minimum calpulse width '''
140 ''' Sets the maximum calpulse amplitude '''
145 ''' Sets the minimum calpulse amplitude '''
150 ''' Sets the flag to ingore the hits without calibration '''
155 ''' Event processor: fill histograms '''
159 for ipulse1, digit1
in enumerate(digits):
162 if (digit1.getHitQuality() != 0
and
168 slotID = digit1.getModuleID()
169 hwchan = digit1.getChannel()
170 for ipulse2, digit2
in enumerate(digits, start=ipulse1 + 1):
174 if (digit2.getHitQuality() != 0
and
179 slotID == digit2.getModuleID()
and
180 hwchan == digit2.getChannel()):
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')
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()
227 ''' Write histograms to file, fills and fits the resolution plots'''
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())
242 projection.Fit(
"gaussFit")
246 tfile = TFile(self.
outname,
'recreate')
254 self.
h_dVdtRising_1.GetXaxis().SetTitle(
"dV/dt [ADC counts / sample]")
256 self.
h_dVdtRising_2.GetXaxis().SetTitle(
"dV/dt [ADC counts / sample]")
258 self.
h_dVdtFalling_1.GetXaxis().SetTitle(
"dV/dt [ADC counts / sample]")
260 self.
h_dVdtFalling_2.GetXaxis().SetTitle(
"dV/dt [ADC counts / sample]")
277 self.
h_DeltaT_RR.GetXaxis().SetTitle(
"#Delta t [ns]")
279 self.
h_DeltaT_FF.GetXaxis().SetTitle(
"#Delta t [ns]")
281 self.
h_DeltaT_FR.GetXaxis().SetTitle(
"#Delta t [ns]")
283 self.
h_DeltaT_RF.GetXaxis().SetTitle(
"#Delta t [ns]")
286 self.
h_DeltaTVSChannel_RR.GetXaxis().SetTitle(
"Global channel number [hwChannel + 512*(slotID-1)]")
289 self.
h_DeltaTVSChannel_FF.GetXaxis().SetTitle(
"Global channel number [hwChannel + 512*(slotID-1)]")
292 self.
h_DeltaTVSChannel_FR.GetXaxis().SetTitle(
"Global channel number [hwChannel + 512*(slotID-1)]")
295 self.
h_DeltaTVSChannel_RF.GetXaxis().SetTitle(
"Global channel number [hwChannel + 512*(slotID-1)]")
299 self.
h_DeltaTVSdVdt_RR.GetXaxis().SetTitle(
"Average of dV/dt on first and second pulse [ACD counts / sample]")
303 self.
h_DeltaTVSdVdt_FF.GetXaxis().SetTitle(
"Average of dV/dt on first and second pulse [ACD counts / sample]")
313print(
'usage: basf2', argvs[0],
'runNumber outfileName')
319runnumbers = sys.argv[5:]
322for runnumber
in runnumbers:
323 files += [f
for f
in glob.glob(str(folder) +
'/*' + str(runnumber) +
'*.sroot')]
327 for runnumber
in runnumbers:
328 files += [f
for f
in glob.glob(str(folder) +
'/*' + str(runnumber) +
'*.root')]
331 print(
"file: " + fname)
333if dbaddress !=
'none':
334 print(
"using local DB " + dbaddress)
335 b2.conditions.append_testing_payloads(dbaddress +
"/localDB.txt")
337 print(
"database not set. Continuing without calibrations")
340b2.set_log_level(b2.LogLevel.ERROR)
343b2.conditions.override_globaltags()
344b2.conditions.append_globaltag(
'online')
347main = b2.create_path()
351 roinput = b2.register_module(
'SeqRootInput')
353 roinput = b2.register_module(
'RootInput')
354roinput.param(
'inputFileNames', files)
355main.add_module(roinput)
358if datatype ==
'pocket':
359 print(
'pocket DAQ data assumed')
360 converter = b2.register_module(
'Convert2RawDet')
361 main.add_module(converter)
364main.add_module(
'TOPGeometryParInitializer')
367unpack = b2.register_module(
'TOPUnpacker')
368main.add_module(unpack)
371featureExtractor = b2.register_module(
'TOPWaveformFeatureExtractor')
372main.add_module(featureExtractor)
375converter = b2.register_module(
'TOPRawDigitConverter')
376if dbaddress ==
'none':
377 print(
"Not using TBC")
378 converter.param(
'useSampleTimeCalibration',
False)
381 converter.param(
'useSampleTimeCalibration',
True)
382converter.param(
'useAsicShiftCalibration',
False)
383converter.param(
'useChannelT0Calibration',
False)
384converter.param(
'useModuleT0Calibration',
False)
385converter.param(
'useCommonT0Calibration',
False)
386converter.param(
'calibrationChannel', -1)
387converter.param(
'lookBackWindows', 28)
388main.add_module(converter)
392resolutionModule.setOutputName(outfile)
393resolutionModule.setMinWidth(0.5)
394resolutionModule.setMaxWidth(3.5)
395resolutionModule.setMinAmp(150)
396resolutionModule.setMaxAmp(999)
397if dbaddress ==
'none':
398 resolutionModule.ignoreNotCalibrated(
False)
400 resolutionModule.ignoreNotCalibrated(
True)
401main.add_module(resolutionModule)
404progress = b2.register_module(
'Progress')
405main.add_module(progress)
412print(b2.statistics(b2.statistics.TERM))
A (simplified) python wrapper for StoreArray.
h_DeltaT_RR
DeltaT rising-rising.
bool m_ignoreNotCalibrated
ignores the hits wthout calibration
h_DeltaTVSChannel_RF
DeltaT rising-falling VS channel.
int m_calpulseMinAmp
minimum amplitude to flag a calpulse candidate
h_WidthVSAmplitude_1
Width VS amplitude, first calibration pulse.
h_DeltaTVSdVdt_RR
DeltaT rising-rising VS average of dV/dt on the first and second pulse.
h_DeltaT_FF
DeltaT falling-falling.
h_DeltaT_RF
DeltaT rising-falling.
h_ResolutionVSdVdt_RR
DeltaT resolution VS average of dV/dt (rising-rising)
h_dVdtFallingDifference
Difference between the dV/dt of the first and the second calpulse, using the rising edges.
h_dVdtFalling_1
dV/dt on the falling edge, first calibration pulse
h_dVdtFalling_2
dV/dt on the falling edge, second calibration pulse
h_dVdtRisingVSdVdtFalling_2
dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse
h_DeltaTVSChannel_RR
DeltaT rising-rising VS channel.
int m_calpulseMaxWidth
maximum width to flag a calpulse candidate
h_dVdtRising_2
dV/dt on the rising edge, second calibration pulse
ignoreNotCalibrated(self, ignoreNotCal)
setMinWidth(self, minWidth)
int m_calpulseMaxAmp
minimum amplitude to flag a calpulse candidate
h_dVdtRisingDifference
Difference between the dV/dt of the first and the second calpulse, using the rising edges.
h_DeltaT_FR
DeltaT falling-rising.
str outname
output root file
h_DeltaTVSdVdt_FF
DeltaT falling-falling VS average of dV/dt on the first and second pulse.
h_dVdtRising_1
dV/dt on the rising edge, first calibration pulse
h_ResolutionVSdVdt_FF
DeltaT resolution VS average of dV/dt (falling-falling)
setOutputName(self, outputname)
float m_calpulseMinWidth
minimum width to flag a calpulse candidate
h_DeltaTVSChannel_FF
DeltaT falling-falling VS channel.
h_dVdtRisingVSdVdtFalling_1
dV/dt on the rising edge VS dV/dt on the falling edge, first calibration pulse
h_DeltaTVSChannel_FR
DeltaT falling-rising VS channel.
setMaxWidth(self, maxWidth)
h_WidthVSAmplitude_2
Width VS amplitude, second calibration pulse.