SND@LHC Software
Loading...
Searching...
No Matches
ScifiCTR.py
Go to the documentation of this file.
1import os,atexit
2import shipunit as u
3import rootUtils as ut
4import ROOT
5import ctypes
6import pickle
7import SndlhcGeo
8from array import array
9def pyExit():
10 print("Make suicide until solution found for freezing")
11 os.system('kill '+str(os.getpid()))
12atexit.register(pyExit)
13
14# supress on-screen histogram drawing
15ROOT.gROOT.SetBatch(ROOT.kTRUE)
16
17A,B=ROOT.TVector3(),ROOT.TVector3()
18
19def FCN(npar, gin, f, par, iflag):
20#calculate chisquare
21 chisq = 0
22 X = ROOT.gROOT.FindObjectAny('commonBlock')
23 for matH in range(nMats):
24 for matV in range(nMats):
25 tdiff = X.GetBinContent(matH*10+matV)
26 d = tdiff - (par[matH] - par[matV+nMats])
27 chisq += d**2
28 f.value = chisq
29 return
30
31def FCNS(npar, gin, f, par, iflag):
32#calculate chisquare
33 chisq = 0
34 X = ROOT.gROOT.FindObjectAny('commonBlock')
35 for s1 in range(1, nStations+1):
36 for s2 in range(s1+1, nStations+1):
37 tdiff = X.GetBinContent(s1*10+s2)
38 d = tdiff - (par[s2-1] - par[s1-1])
39 chisq += d**2
40 f.value = chisq
41 return
42
43class Scifi_CTR(ROOT.FairTask):
44 " time resolution comparing X/Y in same station"
45 def Init(self,monitor):
46 self.M = monitor
47 h = self.M.h
48 self.projs = {1:'V',0:'H'}
49
50 # setup geometry
51 if (options.geoFile).find('../')<0: self.snd_geo = SndlhcGeo.GeoInterface(options.path+options.geoFile)
52 else: self.snd_geo = SndlhcGeo.GeoInterface(options.geoFile[3:])
53 self.MuFilter = self.snd_geo.modules['MuFilter']
54 self.Scifi = self.snd_geo.modules['Scifi']
55
56 self.tag = monitor.iteration
57 tag = self.tag
58 for s in range(1, nStations+1):
59 ut.bookHist(h,'CTR_Scifi'+str(s)+tag,'CTR '+str(s)+tag+'; dt [ns]; ',100,-5.,5.)
60 ut.bookHist(h,'CTR_Scifi_beam'+str(s)+tag,'CTR beam'+str(s)+tag+'; dt [ns]; ',100,-5.,5.)
61 for matH in range(nMats):
62 for matV in range(nMats):
63 ut.bookHist(h,'CTR_Scifi'+str(s*100+10*matH+matV)+tag,'CTR '+str(s)+tag+'; dt [ns]; ',100,-5.,5.)
64 ut.bookHist(h,'CTR_Scifi_beam'+str(s*100+10*matH+matV)+tag,'CTR beam'+str(s)+tag+'; dt [ns]; ',100,-5.,5.)
65
66 for p in self.projs:
67 ut.bookHist(h,'res'+str(s)+self.projs[p],'d; [cm]; ',100,-1.,1.)
68 ut.bookHist(h,'resX'+str(s)+self.projs[p],'d; [cm]; ',100,-1.,1.)
69 ut.bookHist(h,'extrap'+str(s)+self.projs[p],'SiPM distance; L [cm]; ',100,0.,50.)
70
71 for s1 in range(1, nStations):
72 for s2 in range(s1+1, nStations+1):
73 ut.bookHist(h,'CTR_ScifiStation'+str(s1*10+s2)+tag,'CTR station'+str(s1*10+s2)+'; dt [ns]; ',100,-5.,5.)
74 ut.bookHist(h,'CTR_ScifiStation_beam'+str(s1*10+s2)+tag,'CTR station beam'+str(s1*10+s2)+'; dt [ns]; ',100,-5.,5.)
75
76 if self.tag == "v0":
78 for s in range(1, nStations+1):
79 if nMats == 1 :
80 self.tdcScifiStationCalib[s] = [0,{'H':{0:0},'V':{0:0}}]
81 if nMats == 3:
82 self.tdcScifiStationCalib[s] = [0,{'H':{0:0,1:0,2:0},'V':{0:0,1:0,2:0}}]
83
84 else:
85 with open('ScifiTimeAlignment_'+self.tag, 'rb') as fh:
86 self.tdcScifiStationCalib = pickle.load(fh)
87
88 self.V = M.Scifi.GetConfParF("Scifi/signalSpeed")
89
90 def ExecuteEvent(self,event):
91 h = self.M.h
92 tag = self.tag
93 DetID2Key={}
94 for nHit in range(event.Digi_ScifiHits.GetEntries()):
95 DetID2Key[event.Digi_ScifiHits[nHit].GetDetectorID()] = nHit
96
97 for aTrack in self.M.Reco_MuonTracks:
98 if not aTrack.GetUniqueID()==1: continue
99 fitStatus = aTrack.getFitStatus()
100 if not fitStatus.isFitConverged(): continue
101 state = aTrack.getFittedState()
102 mom = state.getMom()
103 slopeX = mom.X()/mom.Z()
104 slopeY = mom.Y()/mom.Z()
105 sortedClusters={}
106#
107 pos = {}
108 for s in range(1, nStations+1): sortedClusters[s] = {'H':[],'V':[]}
109 for nM in range(aTrack.getNumPointsWithMeasurement()):
110 state = aTrack.getFittedState(nM)
111 Meas = aTrack.getPointWithMeasurement(nM)
112 W = Meas.getRawMeasurement()
113 clkey = W.getHitId()
114 aCl = self.M.trackTask.clusScifi[clkey]
115 aHit = event.Digi_ScifiHits[DetID2Key[aCl.GetFirst()]]
116 s = aCl.GetFirst()//1000000
117 aCl.GetPosition(A,B)
118 mat = (aCl.GetFirst()//10000)%10
119 if aHit.isVertical():
120 L = B[1]-state.getPos()[1]
121 sortedClusters[s]['V'].append( [clkey,L,B[0],state.getPos()[1],mat,(A[2]+B[2])/2.] )
122 rc = h['res'+str(s)+'V'].Fill( (A[0]+B[0])/2.-state.getPos()[0])
123 else:
124 L = A[0]-state.getPos()[0]
125 sortedClusters[s]['H'].append( [clkey,L,A[1],state.getPos()[0],mat,(A[2]+B[2])/2.] )
126 rc = h['res'+str(s)+'H'].Fill( (A[1]+B[1])/2.-state.getPos()[1])
127 # find station with exactly 1 x and 1 y cluster:
128 for s in range(1, nStations+1):
129 if not (len(sortedClusters[s]['V']) * len(sortedClusters[s]['H']) ) ==1: continue
130 timeCorr = {}
131 for proj in ['V','H']:
132 clkey = sortedClusters[s][proj][0][0]
133 aCl = self.M.trackTask.clusScifi[clkey]
134 L = sortedClusters[s][proj][0][1]
135 rc = h['extrap'+str(s)+proj].Fill(abs(L))
136 time = aCl.GetTime() # Get time in ns, use fastest TDC of cluster
137 if M.options.check: time = self.M.Scifi.GetCorrectedTime(aCl.GetFirst(),aCl.GetTime(),0)
138 mat = sortedClusters[s][proj][0][4]
139 time-= self.tdcScifiStationCalib[s][1][proj][mat]
140 time-= self.tdcScifiStationCalib[s][0]
141 timeCorr[proj] = time - abs(L)/self.V
142 # print(s,proj,time,L,L/self.V)
143 dt = timeCorr['H'] - timeCorr['V']
144 rc = h['CTR_Scifi'+str(s)+tag].Fill(dt)
145 matH,matV = sortedClusters[s]['H'][0][4],sortedClusters[s]['V'][0][4]
146 rc = h['CTR_Scifi'+str(100*s+10*matH+matV)+tag].Fill(dt)
147 if abs(slopeX)<0.1 and abs(slopeY)<0.1:
148 rc = h['CTR_Scifi_beam'+str(s)+tag].Fill(dt)
149 rc = h['CTR_Scifi_beam'+str(100*s+10*matH+matV)+tag].Fill(dt)
150 dR = sortedClusters[s]['V'][0][2] - sortedClusters[s]['H'][0][3]
151 rc = h['resX'+str(s)+'V'].Fill(dR)
152 dR = sortedClusters[s]['H'][0][2]-sortedClusters[s]['V'][0][3]
153 rc = h['resX'+str(s)+'H'].Fill(dR)
154#
155 if tag!='v0':
156 stationTimes = {}
157 for s in range(1, nStations+1):
158 if not (len(sortedClusters[s]['V']) * len(sortedClusters[s]['H']) ) ==1: continue
159 sTime = 0
160 for proj in ['V','H']:
161 clkey = sortedClusters[s][proj][0][0]
162 aCl = self.M.trackTask.clusScifi[clkey]
163 L = sortedClusters[s][proj][0][1]
164 time = aCl.GetTime() # Get time in ns, use fastest TDC of cluster
165 if M.options.check: time = self.M.Scifi.GetCorrectedTime(aCl.GetFirst(),aCl.GetTime(),0)
166 mat = sortedClusters[s][proj][0][4]
167 time-= self.tdcScifiStationCalib[s][1][proj][mat]
168 time-= self.tdcScifiStationCalib[s][0]
169 time-= abs(L)/self.V
170 sTime += time
171 stationTimes[s] = [sTime/2.,(sortedClusters[s]['H'][0][5] + sortedClusters[s]['V'][0][5])/2.]
172 for s1 in range(1, nStations):
173 if not s1 in stationTimes: continue
174 for s2 in range(s1+1, nStations+1):
175 if not s2 in stationTimes: continue
176 dT = stationTimes[s2][0] - stationTimes[s1][0]
177# correct for slope
178 dZ = stationTimes[s2][1] - stationTimes[s1][1]
179 dL = dZ * ROOT.TMath.Sqrt( slopeX**2+slopeY**2+1 )
180 if slopeY>0.1: dL = -dL # cosmics from the back
181 dT -= dL / u.speedOfLight
182 rc = h['CTR_ScifiStation'+str(s1*10+s2)+tag].Fill(dT)
183 if abs(slopeX)<0.1 and abs(slopeY)<0.1:
184 rc = h['CTR_ScifiStation_beam'+str(s1*10+s2)+tag].Fill(dT)
185
186 def Plot(self):
187 h = self.M.h
188 tag = self.tag
189 for b in ['','_beam']:
190 ut.bookCanvas(h,'CTR'+b+tag,'CTR'+tag,1800,1200,3,2)
191 for s in range(1, nStations+1):
192 tc = h['CTR'+b+tag].cd(s)
193 histo = h['CTR_Scifi'+b+str(s)+tag]
194 rc = histo.Fit('gaus','SQ')
195 tc.Update()
196 stats = histo.FindObject('stats')
197 stats.SetOptFit(1111111)
198 stats.SetX1NDC(0.62)
199 stats.SetY1NDC(0.63)
200 stats.SetX2NDC(0.98)
201 stats.SetY2NDC(0.94)
202 histo.Draw()
203
204 for s in range(1, nStations+1):
205 ut.bookCanvas(h,'CTRM'+str(s)+b+tag,'CTR per mat combi',1800,1200,3,3)
206 j = 1
207 for matH in range(nMats):
208 for matV in range(nMats):
209 tc = h['CTRM'+str(s)+b+tag].cd(j)
210 j+=1
211 histo = h['CTR_Scifi'+b+str(100*s+10*matH+matV)+tag]
212 histo.Fit('gaus','SQ')
213 tc.Update()
214 stats = histo.FindObject('stats')
215 stats.SetOptFit(1111111)
216 stats.SetX1NDC(0.62)
217 stats.SetY1NDC(0.63)
218 stats.SetX2NDC(0.98)
219 stats.SetY2NDC(0.94)
220 histo.Draw()
221
222 if tag!='v0':
223 for b in ['','_beam']:
224 ut.bookCanvas(h,'CTRS'+b,'CTR per station combination'+tag,1800,1200,5,2)
225 k=1
226 for s1 in range(1, nStations):
227 for s2 in range(s1+1, nStations+1):
228 tc = h['CTRS'+b].cd(k)
229 k+=1
230 histo = h['CTR_ScifiStation'+b+str(s1*10+s2)+tag]
231 histo.Draw()
232 histo.Fit('gaus','SQ')
233 tc.Update()
234 stats = histo.FindObject('stats')
235 stats.SetOptFit(1111111)
236 stats.SetX1NDC(0.62)
237 stats.SetY1NDC(0.63)
238 stats.SetX2NDC(0.98)
239 stats.SetY2NDC(0.94)
240 histo.Draw()
241
243 self.meanAndSigma = {}
244 h = self.M.h
245 tag = self.tag
246 ut.bookHist(h,'cor')
247 for b in ['','_beam']:
248 self.meanAndSigma[b]={}
249 for s in range(1, nStations+1):
250 for matH in range(nMats):
251 for matV in range(nMats):
252 key = 100*s+10*matH+matV
253 histo = h['CTR_Scifi'+b+str(key)+tag]
254 Fun = histo.GetFunction('gaus')
255 self.meanAndSigma[b][key] = [Fun.GetParameter(1),Fun.GetParameter(2)]
256
257 def minimize(self,b=""):
258 h = self.M.h
260 ut.bookHist(h,'commonBlock','',100,0.,100.)
261
262 for s in range(1, nStations+1):
263 for matH in range(nMats):
264 for matV in range(nMats):
265 key = 100*s+10*matH+matV
266 dt = self.meanAndSigma[b][key][0]
267 h['commonBlock'].SetBinContent(matH*10+matV,dt)
268 npar = 2*nMats
269 ierflg = ctypes.c_int(0)
270 vstart = array('d',[0]*npar)
271 gMinuit = ROOT.TMinuit(npar)
272 gMinuit.SetFCN(FCN)
273 gMinuit.SetErrorDef(1.0)
274 gMinuit.SetMaxIterations(10000)
275 err = 1E-3
276 p = 0
277 for proj in ['H','V']:
278 for m in range(nMats):
279 name = "s"+str(s)+proj+str(m)
280 gMinuit.mnparm(p, name, vstart[p], err, 0.,0.,ierflg)
281 p+=1
282 gMinuit.FixParameter(0)
283 strat = array('d',[0])
284 gMinuit.mnexcm("SET STR",strat,1,ierflg) # 0 faster, 2 more reliable
285 gMinuit.mnexcm("SIMPLEX",vstart,npar,ierflg)
286 gMinuit.mnexcm("MIGRAD",vstart,npar,ierflg)
287
288 cor = ctypes.c_double(0)
289 ecor = ctypes.c_double(0)
290 p = 0
291 for proj in ['H','V']:
292 for m in range(nMats):
293 rc = gMinuit.GetParameter(p,cor,ecor)
294 p+=1
295 self.tdcScifiStationCalib[s][1][proj][m] = cor.value
296
297 with open('ScifiTimeAlignment_v1', 'wb') as fh:
298 pickle.dump(self.tdcScifiStationCalib, fh)
299
300 def minimizeStation(self,b=""):
301 h = self.M.h
302 tag = self.tag
304 self.meanAndSigmaStation[b]={}
305 for s1 in range(1, nStations+1):
306 for s2 in range(s1+1, nStations+1):
307 key = s1*10+s2
308 histo = h[ 'CTR_ScifiStation'+b+str(key)+tag]
309 Fun = histo.GetFunction('gaus')
310 self.meanAndSigmaStation[b][key] = [Fun.GetParameter(1),Fun.GetParameter(2)]
311 for s1 in range(1, nStations+1):
312 for s2 in range(s1+1, nStations+1):
313 key = 10*s1+s2
314 dt = self.meanAndSigmaStation[b][key][0]
315 h['commonBlock'].SetBinContent(key,dt)
316 npar = 5
317 ierflg = ctypes.c_int(0)
318 vstart = array('d',[0,-0.3365,-0.975,0.11,0.239])
319 gMinuit = ROOT.TMinuit(npar)
320 gMinuit.SetFCN(FCNS)
321 gMinuit.SetErrorDef(1.0)
322 gMinuit.SetMaxIterations(10000)
323 err = 1E-3
324 p = 0
325 for s in range(1, nStations+1):
326 name = "station"+str(s)
327 gMinuit.mnparm(p, name, vstart[p], err, 0.,0.,ierflg)
328 p+=1
329 gMinuit.FixParameter(0)
330 strat = array('d',[0])
331 gMinuit.mnexcm("SET STR",strat,1,ierflg) # 0 faster, 2 more reliable
332 gMinuit.mnexcm("SIMPLEX",vstart,npar,ierflg)
333 gMinuit.mnexcm("MIGRAD",vstart,npar,ierflg)
334
335 cor = ctypes.c_double(0)
336 ecor = ctypes.c_double(0)
337 for s in range(1, nStations+1):
338 rc = gMinuit.GetParameter(s-1,cor,ecor)
339 self.tdcScifiStationCalib[s][0] = cor.value
340
341 with open('ScifiTimeAlignment_v2', 'wb') as fh:
342 pickle.dump(self.tdcScifiStationCalib, fh)
343
344class Scifi_TimeOfTracks(ROOT.FairTask):
345 " dt versus dz"
346 def Init(self,monitor):
347 self.M = monitor
348 self.trackTask = self.M.FairTasks['simpleTracking']
349 h = self.M.h
350 self.V = M.Scifi.GetConfParF("Scifi/signalSpeed")
351 ut.bookHist(M.h,'dTvsZ','dT versus dL; dL/0.5cm [cm];dT/100ps [ns]',140,0.,70.,120,-6.,6.)
352 ut.bookHist(M.h,'dTvsZ_beam','dT versus dL, beam dL/0.5cm [cm];dT/100ps [ns]',140,0.,70.,120,-6.,6.)
353 with open('ScifiTimeAlignment_v2', 'rb') as fh:
354 self.tdcScifiStationCalib = pickle.load(fh)
355
356 def ExecuteEvent(self,event):
357 h = self.M.h
358 for aTrack in self.M.Reco_MuonTracks:
359 if not aTrack.GetUniqueID()==1: continue
360 fitStatus = aTrack.getFitStatus()
361 if not fitStatus.isFitConverged(): continue
362 state = aTrack.getFittedState()
363 mom = state.getMom()
364 slopeX = mom.X()/mom.Z()
365 slopeY = mom.Y()/mom.Z()
366 sortedClusters={}
367 DetID2Key={}
368 for nHit in range(event.Digi_ScifiHits.GetEntries()):
369 DetID2Key[event.Digi_ScifiHits[nHit].GetDetectorID()] = nHit
370#
371 pos = {}
372 for s in range(1, nStations+1): sortedClusters[s] = {'H':[],'V':[]}
373 for nM in range(aTrack.getNumPointsWithMeasurement()):
374 state = aTrack.getFittedState(nM)
375 Meas = aTrack.getPointWithMeasurement(nM)
376 W = Meas.getRawMeasurement()
377 clkey = W.getHitId()
378 aCl = self.trackTask.clusScifi[clkey]
379 aHit = event.Digi_ScifiHits[DetID2Key[aCl.GetFirst()]]
380 s = aCl.GetFirst()//1000000
381 aCl.GetPosition(A,B)
382 mat = (aCl.GetFirst()//10000)%10
383
384 if aHit.isVertical():
385 L = B[1]-state.getPos()[1]
386 sortedClusters[s]['V'].append( [clkey,L,B[0],state.getPos()[1],mat,(A[2]+B[2])/2.] )
387 else:
388 L = A[0]-state.getPos()[0]
389 sortedClusters[s]['H'].append( [clkey,L,A[1],state.getPos()[0],mat,(A[2]+B[2])/2.] )
390 stationTimes = {}
391 for s in range(1, nStations+1):
392 if not (len(sortedClusters[s]['V']) * len(sortedClusters[s]['H']) ) ==1: continue
393 sTime = 0
394 for proj in ['V','H']:
395 clkey = sortedClusters[s][proj][0][0]
396 aCl = self.trackTask.clusScifi[clkey]
397 L = sortedClusters[s][proj][0][1]
398 time = aCl.GetTime() # Get time in ns, use fastest TDC of cluster
399 mat = sortedClusters[s][proj][0][4]
400 time-= self.tdcScifiStationCalib[s][1][proj][mat] # correct as function of station / projection / mat
401 time-= self.tdcScifiStationCalib[s][0] # internal station calibration
402 time-= abs(L)/self.V
403
404# cross check
405 # corTime = self.M.Scifi.GetCorrectedTime(aCl.GetFirst(),aCl.GetTime(),abs(L) )
406 # print('test',aCl.GetTime(),time,corTime)
407
408 sTime += time
409 stationTimes[s] = [sTime/2.,(sortedClusters[s]['H'][0][5] + sortedClusters[s]['V'][0][5])/2.]
410 # require station 1 to be present which defines T0
411 if not 1 in stationTimes: continue
412
413 T0 = stationTimes[1][0]
414 Z0 = stationTimes[1][1]
415 for s in stationTimes:
416 if s == 1: continue
417 dZ = stationTimes[s][1]-Z0
418 dT = stationTimes[s][0]-T0
419 dL = dZ * ROOT.TMath.Sqrt( slopeX**2+slopeY**2+1 )
420 # if slopeY>0.1: dL = -dL # cosmics from the back
421 # dT -= dL / u.speedOfLight
422 rc = h['dTvsZ'].Fill(dL,dT)
423 if abs(slopeX)<0.1 and abs(slopeY)<0.1:
424 rc = h['dTvsZ_beam'].Fill(dL,dT)
425
426
427
428 def Plot(self):
429 h = self.M.h
430 ut.bookCanvas(h,'Tscifi','',1800,1200,2,1)
431 h['Tscifi'].cd(1)
432 h['dTvsZ'].Draw('colz')
433 h['Tscifi'].cd(2)
434 h['dTvsZ_beam'].Draw('colz')
435
436class histStore():
437 h = {}
438 xCheck = ''
439
440if __name__ == '__main__':
441 import Monitor
442 import SndlhcTracking
443 from argparse import ArgumentParser
444 parser = ArgumentParser()
445 parser.add_argument("--server", dest="server", help="xrootd server",default=os.environ["EOSSHIP"])
446 parser.add_argument("-r", "--runNumber", dest="runNumber", help="run number", type=int,default=-1)
447 parser.add_argument("-p", "--path", dest="path", help="path to data",required=False,default="/eos/experiment/sndlhc/convertedData/physics/2022/")
448 parser.add_argument("-P", "--partition", dest="partition", help="partition of data", type=int,required=False,default=-1)
449 parser.add_argument("-f", "--inputFile", dest="fname", help="file name", type=str,default=None,required=False)
450 parser.add_argument("-g", "--geoFile", dest="geoFile", help="file name", type=str,default="geofile_sndlhc_TI18_V0_2022.root",required=False)
451 parser.add_argument("-c", "--command", dest="command", help="command",required=False,default="")
452 parser.add_argument("-o", dest="check", help="use corrected times",action='store_true',default=False)
453 parser.add_argument("-M", "--online", dest="online", help="online mode",default=False,action='store_true')
454 parser.add_argument("--interactive", dest="interactive", action='store_true',default=False)
455 parser.add_argument("-n", "--nEvents", dest="nEvents", help="number of events", default=-1,type=int)
456
457 options = parser.parse_args()
458 options.trackType = 'Scifi'
459
460 FairTasks = []
462 trackTask.SetName('simpleTracking')
463 FairTasks.append(trackTask)
464 M = Monitor.Monitoring(options,FairTasks)
465 nMats = M.Scifi.GetConfParI("Scifi/nmats")
466 nStations = M.Scifi.GetConfParI("Scifi/nscifi")
467
468 if options.nEvents <0: options.nEvents = M.GetEntries()
469 if options.command == "full":
470 task = Scifi_CTR()
471 M.iteration = 'v0'
472 task.Init(M)
473 for n in range(options.nEvents):
474 event = M.GetEvent(n)
475 task.ExecuteEvent(event)
476 task.Plot()
477 task.minimize(b="") # all tracks, not yet enough stats for beam tracks
478# cross check and station alignment
479 M.iteration = 'v1'
480 task.Init(M)
481 for n in range(options.nEvents):
482 event = M.GetEvent(n)
483 task.ExecuteEvent(event)
484 task.Plot()
485# nw minimize the station alignments
486 task.minimizeStation(b="")
487
488# cross check
489 M.iteration = 'v2'
490 task.Init(M)
491 for n in range(options.nEvents):
492 event = M.GetEvent(n)
493 task.ExecuteEvent(event)
494 task.Plot()
495 ut.writeHists(M.h,'ScifiTimeCalibration.root',plusCanvas=True)
496
497 for s in range(1, nStations+1):
498 C = task.tdcScifiStationCalib[s]
499 if nMats == 1 :
500 print ("station %1i %5.2F H: %5.2F V: %5.2F "%( s, C[0],C[1]['H'][0],C[1]['V'][0] ) )
501 if nMats == 3 :
502 print ("station %1i %5.2F H: %5.2F %5.2F %5.2F V: %5.2F %5.2F %5.2F "%( s, C[0],C[1]['H'][0], C[1]['H'][1], C[1]['H'][2],C[1]['V'][0], C[1]['V'][1], C[1]['V'][2] ) )
503 for s in range(1, nStations+1):
504 C = task.tdcScifiStationCalib[s]
505 if nMats == 1 :
506 print ("station %1i %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns "%( s, C[0],C[1]['H'][0],C[1]['V'][0] ) )
507 if nMats == 3 :
508 print ("station %1i %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns "%( s, C[0],C[1]['H'][0], C[1]['H'][1], C[1]['H'][2],C[1]['V'][0], C[1]['V'][1], C[1]['V'][2] ) )
509
510if options.command == "full" or options.command == "check":
511# final test
513 taskT.Init(M)
514 for n in range(options.nEvents):
515 event = M.GetEvent(n)
516 taskT.ExecuteEvent(event)
517 taskT.Plot()
518 for s in range(1, nStations+1):
519 C = taskT.tdcScifiStationCalib[s]
520 if nMats == 1 :
521 print ("station %1i %5.2F H: %5.2F V: %5.2F "%( s, C[0],C[1]['H'][0],C[1]['V'][0] ) )
522 if nMats == 3 :
523 print ("station %1i %5.2F H: %5.2F %5.2F %5.2F V: %5.2F %5.2F %5.2F "%( s, C[0],C[1]['H'][0], C[1]['H'][1], C[1]['H'][2],C[1]['V'][0], C[1]['V'][1], C[1]['V'][2] ) )
524 for s in range(1, nStations+1):
525 C = taskT.tdcScifiStationCalib[s]
526 if nMats == 1 :
527 print ("station %1i %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns "%( s, C[0],C[1]['H'][0],C[1]['V'][0] ) )
528 if nMats == 3 :
529 print ("station %1i %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns, %5.3F*u.ns "%( s, C[0],C[1]['H'][0], C[1]['H'][1], C[1]['H'][2],C[1]['V'][0], C[1]['V'][1], C[1]['V'][2] ) )
530
Init(self, monitor)
Definition ScifiCTR.py:45
ExecuteEvent(self, event)
Definition ScifiCTR.py:90
minimizeStation(self, b="")
Definition ScifiCTR.py:300
ExtractMeanAndSigma(self)
Definition ScifiCTR.py:242
Double_t GetCorrectedTime(Int_t fDetectorID, Double_t rawTime, Double_t L)
Definition Scifi.cxx:517
FCN(npar, gin, f, par, iflag)
Definition ScifiCTR.py:19
FCNS(npar, gin, f, par, iflag)
Definition ScifiCTR.py:31
pyExit()
Definition ScifiCTR.py:9