Belle II Software  release-08-01-10
studyLaserLight.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # ------------------------------------------------------------------------
13 # Module to produce the laser time histograms that serve as input to the laser
14 # resolution fitter
15 #
16 # usage: basf2 studyLaserResolution.py dbaddress (path|none) type (local|pocket|root)
17 # output_name.root path_to_sroot run1 run2 ... runN
18 # The run number accepts wildcards
19 # ------------------------------------------------------------------------
20 
21 
22 import basf2 as b2
23 import sys
24 import glob
25 from ROOT import Belle2
26 from ROOT import TFile, TH2F
27 from laserResolutionTools import fitLaserResolution, plotLaserResolution
28 
29 
30 class TOPLaserHistogrammerModule(b2.Module):
31 
32  ''' Module to study resolution and performance of the top laser calibration.'''
33 
34 
35  h_LaserTimingVSChannel = TH2F(
36  'LaserTimingVSChannel',
37  'Laser timing in as function of the channel number',
38  512 * 16,
39  0,
40  512 * 16,
41  10000,
42  0.,
43  100) # 10 ps binning
44 
45  h_LaserTimingVSChannelOneSlot = TH2F(
46  'LaserTimingVSChannelOneSlot',
47  'Laser timing in as function of the channel number',
48  512,
49  0,
50  512,
51  10000,
52  0.,
53  100) # 10 ps binning
54 
55 
56  h_crossOccupancy = [[TH2F(
57  'crossOccupancy_' + str(slotA) + '_' + str(slotB),
58  ' ',
59  200,
60  0,
61  200,
62  200,
63  0.,
64  200) for slotA in range(16)] for slotB in range(16)]
65 
66 
67  outname = 'outStudyLaserResolution.root'
68 
69 
70  m_ignoreNotCalibrated = True
71 
72  m_runOnData = True
73 
74 
75  m_maxWidth = 3.
76 
77  m_minWidth = 0.5
78 
79  m_maxAmp = 700.
80 
81  m_minAmp = 250.
82 
83  m_mcCorrectionsFile = '/group/belle2/group/detector/TOP/calibration/MCreferences/t0MC.root'
84 
85  m_MCPeaks = [[]]
86 
87  def setOutputName(self, outputname):
88  ''' Sets the output file name '''
89 
90  self.outnameoutnameoutname = outputname
91 
92  def setMaxWidth(self, maxWidth):
93  ''' Sets the maximum calpulse width '''
94 
95  self.m_maxWidthm_maxWidthm_maxWidth = maxWidth
96 
97  def setMinWidth(self, minWidth):
98  ''' Sets the minimum calpulse width '''
99 
100  self.m_minWidthm_minWidthm_minWidth = minWidth
101 
102  def setMaxAmp(self, maxAmp):
103  ''' Sets the maximum calpulse amplitude '''
104 
105  self.m_maxAmpm_maxAmpm_maxAmp = maxAmp
106 
107  def setMinAmp(self, minAmp):
108  ''' Sets the minimum calpulse amplitude '''
109 
110  self.m_minAmpm_minAmpm_minAmp = minAmp
111 
112  def setMCCorrectionsFile(self, MCfile):
113  ''' Sets the file containing the MC correction'''
114 
115  self.m_mcCorrectionsFilem_mcCorrectionsFilem_mcCorrectionsFile = MCfile
116 
117  def ignoreNotCalibrated(self, ignoreNotCal):
118  ''' Sets the flag to ingore the hits without calibration '''
119 
120  self.m_ignoreNotCalibratedm_ignoreNotCalibratedm_ignoreNotCalibrated = ignoreNotCal
121 
122  def event(self):
123  ''' Event processor: fill histograms '''
124 
125  digits = Belle2.PyStoreArray('TOPDigits')
126  nhits = [0 for i in range(16)]
127  for digit in digits:
128  if(not self.ignoreNotCalibratedignoreNotCalibrated and not digit.isTimeBaseCalibrated()):
129  continue
130  if (digit.getHitQuality() == 1 and
131  digit.getPulseWidth() > self.m_minWidthm_minWidthm_minWidth and digit.getPulseWidth() < self.m_maxWidthm_maxWidthm_maxWidth and
132  digit.getPulseHeight() > self.m_minAmpm_minAmpm_minAmp and digit.getPulseHeight() < self.m_maxAmpm_maxAmpm_maxAmp):
133  slotID = digit.getModuleID()
134  hwchan = digit.getChannel()
135  self.h_LaserTimingVSChannelh_LaserTimingVSChannel.Fill(512 * (slotID - 1) + hwchan, digit.getTime())
136  simhits = digit.getRelationsWith('TOPSimHits')
137  nhits[slotID - 1] = nhits[slotID - 1] + 1
138  for simhit in simhits:
139  self.h_LaserTimingVSChannelOneSloth_LaserTimingVSChannelOneSlot.Fill(hwchan, simhit.getTime())
140  for slotA in range(16):
141  for slotB in range(16):
142  self.h_crossOccupancyh_crossOccupancy[slotA][slotB].Fill(nhits[slotA], nhits[slotB])
143 
144  def terminate(self):
145  ''' Write histograms to file, fills and fits the resolution plots'''
146  tfile = TFile(self.outnameoutnameoutname, 'recreate')
147  self.h_LaserTimingVSChannelh_LaserTimingVSChannel.Write()
148  self.h_LaserTimingVSChannelOneSloth_LaserTimingVSChannelOneSlot.Write()
149  for slotA in range(16):
150  for slotB in range(16):
151  self.h_crossOccupancyh_crossOccupancy[slotA][slotB].Write()
152  tfile.Close()
153 
154 
155 argvs = sys.argv
156 print('usage: basf2', argvs[0], 'runNumber outfileName')
157 
158 lookBack = argvs[1] # lookbackWindows setting
159 dbaddress = argvs[2] # path to the calibration DB (absolute path or 'none')
160 datatype = argvs[3] # data type ('pocket' or 'local' or 'root', if root files have to be analyzed)
161 outfile = argvs[4] # output name
162 folder = argvs[5] # folder containing the sroot files
163 runnumbers = sys.argv[6:] # run numbers
164 
165 files = []
166 for runnumber in runnumbers:
167  if datatype == 'root':
168  files += [f for f in glob.glob(str(folder) + '/*' + str(runnumber) + '*.root')]
169  else:
170  # print('running on sroot files in folder ' + glob.glob(str(folder))
171  files += [f for f in glob.glob(str(folder) + '/*' + str(runnumber) + '*.sroot')]
172 for fname in files:
173  print("file: " + fname)
174 
175 if dbaddress != 'none':
176  print("using local DB " + dbaddress)
177  b2.conditions.append_testing_payloads(dbaddress + "/localDB.txt")
178 else:
179  print("database not set. Continuing without calibrations")
180 
181 # Suppress messages and warnings during processing
182 b2.set_log_level(b2.LogLevel.ERROR)
183 
184 # Define a global tag (note: the one given bellow can be out-dated!)
185 b2.conditions.append_globaltag('data_reprocessing_proc8')
186 
187 # Create path
188 main = b2.create_path()
189 
190 # input
191 if datatype == 'root':
192  roinput = b2.register_module('RootInput')
193  roinput.param('inputFileNames', files)
194  main.add_module(roinput)
195 else:
196  roinput = b2.register_module('SeqRootInput')
197  roinput.param('inputFileNames', files)
198  main.add_module(roinput)
199 
200 # conversion from RawCOPPER or RawDataBlock to RawDetector objects
201 if datatype == 'pocket':
202  print('pocket DAQ data assumed')
203  converter = b2.register_module('Convert2RawDet')
204  main.add_module(converter)
205 
206 # Initialize TOP geometry parameters (creation of Geant geometry is not needed)
207 main.add_module('TOPGeometryParInitializer')
208 
209 if datatype != 'root':
210  # Unpacking (format auto detection works now)
211  unpack = b2.register_module('TOPUnpacker')
212  main.add_module(unpack)
213 
214  # Add multiple hits by running feature extraction offline
215  featureExtractor = b2.register_module('TOPWaveformFeatureExtractor')
216  main.add_module(featureExtractor)
217 
218  # Convert to TOPDigits
219  converter = b2.register_module('TOPRawDigitConverter')
220  if dbaddress == 'none':
221  print("Not using Calibrations")
222  converter.param('useSampleTimeCalibration', False)
223  converter.param('useAsicShiftCalibration', False)
224  converter.param('useChannelT0Calibration', False)
225  converter.param('useModuleT0Calibration', False)
226  else:
227  print("Using Calibrations")
228  converter.param('useSampleTimeCalibration', True)
229  converter.param('useAsicShiftCalibration', True)
230  converter.param('useChannelT0Calibration', True)
231  converter.param('useModuleT0Calibration', True)
232  converter.param('useCommonT0Calibration', False)
233  converter.param('calibrationChannel', -1) # do not specify the calpulse channel
234  converter.param('lookBackWindows', int(lookBack)) # in number of windows
235  main.add_module(converter)
236 
237 # resolution plotter
238 resolutionModule = TOPLaserHistogrammerModule()
239 resolutionModule.setOutputName(outfile)
240 resolutionModule.setMinWidth(0.1) # good TOPDigit selection
241 resolutionModule.setMaxWidth(999) # good TOPDigit selection
242 resolutionModule.setMinAmp(100) # good TOPDigit selection
243 resolutionModule.setMaxAmp(999) # good TOPDigit selection
244 if dbaddress == 'none':
245  resolutionModule.ignoreNotCalibrated(True)
246 else:
247  resolutionModule.ignoreNotCalibrated(False)
248 main.add_module(resolutionModule)
249 
250 # Show progress of processing
251 progress = b2.register_module('Progress')
252 main.add_module(progress)
253 
254 # Process events
255 b2.process(main)
256 
257 # Print call statistics
258 print(b2.statistics)
259 print(b2.statistics(b2.statistics.TERM))
260 
261 
262 fitLaserResolution(dataFile=outfile, outputFile='laserResolutionResults.root', pdfType='cb', saveFits=True)
263 plotLaserResolution()
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
int m_maxWidth
maximum width to accept a TOPDigit
bool m_ignoreNotCalibrated
ignores the hits wthout calibration
string m_mcCorrectionsFile
root file with the MC distribution of the laser light, to get the light path corrections
int m_minAmp
minimum amplitude to accept a TOPDigit
int m_maxAmp
minimum amplitude to accept a TOPDigit
h_LaserTimingVSChannelOneSlot
Laser timing in as function of the channel number.
h_LaserTimingVSChannel
Width as function of the sample number in each channel.
float m_minWidth
minimum width to accept a TOPDigit