GeneralBrokenLines V03-01-02
gblpy
gbltst.py
Go to the documentation of this file.
1'''
2Simple Test Program for General Broken Lines.
3
4Created on Jul 27, 2011
5
6@author: kleinwrt
7'''
8
9
29
30import numpy as np
31import math
32import time
33from gblfit import GblPoint, GblTrajectory
34#
35
36
37
56
57#
58 np.random.seed(47117)
59
60 nTry = 1000 #: number of tries
61 nLayer = 10 #: number of detector layers
62 useThickScatterer = False
63 print " Gbltst $Id$ ", nTry, nLayer, useThickScatterer
64 start = time.clock()
65# track direction
66 sinLambda = 0.3
67 cosLambda = math.sqrt(1.0 - sinLambda ** 2)
68 sinPhi = 0.
69 cosPhi = math.sqrt(1.0 - sinPhi ** 2)
70# Curvilinear system: track direction T, U = Z x T / |Z x T|, V = T x U
71 tuvDir = np.array([[cosLambda * cosPhi, cosLambda * sinPhi, sinLambda], \
72 [-sinPhi, cosPhi, 0.], \
73 [-sinLambda * cosPhi, -sinLambda * sinPhi, cosLambda]])
74# measurement resolution
75 measErr = np.array([ 0.001, 0.001]) # 10 mu
76 measPrec = 1.0 / measErr ** 2
77# scattering error
78 scatErr = 0.001 # 1 mrad
79# scattering cov matrix for thick scatterers
80 scatCov = np.zeros((4, 4))
81# RMS of CurviLinear track parameters (Q/P, slopes, offsets)
82 clErr = np.array([0.001, -0.1, 0.2, -0.15, 0.25])
83 # precision matrix for external seed (in local system)
84 locSeed = None
85 seedLabel = 0 # label of point with seed
86 if seedLabel != 0:
87 print " external seed at label ", seedLabel
88#
89 bfac = 0.2998 # Bz*c for Bz=1
90 step = 1.5 / cosLambda # constant steps in RPhi
91#
92 Chi2Sum = 0.
93 NdfSum = 0
94 LostSum = 0.
95#
96 binaryFile = open("milleBinaryISN.dat", "wb")
97#
98 for iTry in range(nTry):
99# generate (CurviLinear) track parameters
100 clNorm = np.random.normal(0., 1., 5)
101 clPar = clErr * clNorm
102 # covariance matrix
103 clCov = np.eye(5)
104 for i in range(5):
105 clCov[i, i] = clErr[i] ** 2
106# arclength
107 s = 0.
108# point-to-point jacobian (from previous point)
109 jacPointToPoint = np.eye(5)
110# additional (local or global) derivatives
111 addDer = np.array([[1.0], [0.0]])
112 labGlobal = np.array([[4711], [4711]])
113# create trajectory
114 traj = GblTrajectory(bfac != 0.)
115
116# at previous point: transformation from local to curvilinear system
117 oldL2c = np.eye(5)
118
119 for iLayer in range(nLayer):
120 #print " layer ", iLayer
121# measurement directions (J,K) from stereo angle
122 sinStereo = (0. if iLayer % 2 == 0 else 0.1)
123 cosStereo = math.sqrt(1.0 - sinStereo ** 2)
124# (orthogonal) measurement system: I, J ,K
125 ijkDir = np.array([[1., 0., 0.], \
126 [0., cosStereo, sinStereo], \
127 [0., -sinStereo, cosStereo]])
128# local system: measurement or curvilinear
129 #local = gblMeasSystem(ijkDir, tuvDir)
130 local = gblCurviSystem(ijkDir, tuvDir)
131# projections
132 proL2m = local.getTransLocalToMeas()
133 proL2c = local.getTransLocalToCurvi()
134 proC2l = np.linalg.inv(proL2c)
135# projection curvilinear to measurement directions
136 proC2m = local.getTransCurviToMeas()
137# measurement - prediction in measurement system with error
138 measNorm = np.random.normal(0., 1., 2)
139 meas = np.dot(proC2m, clPar[3:5]) + measErr * measNorm
140# jacobian is calculated in curvilinear system and transformed
141 jac = np.dot(proC2l, np.dot(jacPointToPoint, oldL2c))
142# point with (independent) measurements (in measurement system)
143 point = GblPoint(jac)
144 point.addMeasurement([proL2m, meas, measPrec])
145# additional local parameters?
146# point.addLocals(addDer)
147# additional global parameters?
148 point.addGlobals(labGlobal, addDer)
149# add thick scatterer
150 if useThickScatterer:
151 if scatCov[0, 0] > 0:
152 scat = np.zeros(4)
153 scatPrec = np.linalg.inv(scatCov)
154 #print " scatCov ", iLayer, scatCov
155 #print " scatPrec ", iLayer, scatPrec
156 point.addScatterer([scat, scatPrec])
157 scatCov = np.zeros((4, 4))
158# add 'noise', 10%
159 measNorm = np.random.uniform(-5., 5., 2)
160 if np.random.uniform() < 0.:
161 meas = np.dot(proC2m, clPar[3:5]) + measErr * measNorm
162 print " noise ", iLayer
163 point.addMeasurement([proL2m, meas, measPrec])
164 point.addGlobals(labGlobal, addDer)
165# locDer flips sign every measurement
166 addDer = -addDer
167# add point to trajectory
168 iLabel = traj.addPoint(point)
169 #print " meas. ", iLabel, s
170 if iLabel == abs(seedLabel):
171 clSeed = np.linalg.inv(clCov)
172 locSeed = np.dot(proL2c, np.dot(clSeed, proL2c.T))
173# propagate to scatterer
174 jacPointToPoint = gblSimpleJacobian(step, cosLambda, bfac)
175 clPar = np.dot(jacPointToPoint, clPar)
176 clCov = np.dot(jacPointToPoint, np.dot(clCov, jacPointToPoint.T))
177 s += step
178# propagate scatCov
179 scatCov = np.dot (jacPointToPoint[1:, 1:], np.dot(scatCov, jacPointToPoint[1:, 1:].T))
180 if (iLayer < nLayer - 1):
181 scat = np.array([0., 0.])
182# point with scatterer
183 jac = np.dot(proC2l, np.dot(jacPointToPoint, proL2c))
184 point = GblPoint(jac)
185# add thin scatterer
186 if not useThickScatterer:
187 if scatErr > 0:
188 scatP = local.getScatPrecision(scatErr)
189 point.addScatterer([scat, scatP])
190 iLabel = traj.addPoint(point)
191 #print " scat1 ", iLabel, s
192 if iLabel == abs(seedLabel):
193 clSeed = np.linalg.inv(clCov)
194 locSeed = np.dot(proL2c, np.dot(clSeed, proL2c.T))
195
196# scatter a little
197 scatNorm = np.random.normal(0., 1., 2)
198 clPar[1:3] = clPar[1:3] + scatErr * scatNorm
199 scatCov[0, 0] += scatErr ** 2
200 scatCov[1, 1] += scatErr ** 2
201# propagate to 2nd scatterer
202 jacPointToPoint = gblSimpleJacobian(step, cosLambda, bfac)
203 clPar = np.dot(jacPointToPoint, clPar)
204 clCov = np.dot(jacPointToPoint, np.dot(clCov, jacPointToPoint.T))
205 s += step
206# propagate scatCov
207 scatCov = np.dot (jacPointToPoint[1:, 1:], np.dot(scatCov, jacPointToPoint[1:, 1:].T))
208 if (iLayer < nLayer - 1):
209 scat = np.array([0., 0.])
210# point with scatterer
211 jac = np.dot(proC2l, np.dot(jacPointToPoint, proL2c))
212 point = GblPoint(jac)
213# add thin scatterer
214 if not useThickScatterer:
215 if scatErr > 0:
216 scatP = local.getScatPrecision(scatErr)
217 point.addScatterer([scat, scatP])
218 iLabel = traj.addPoint(point)
219 #print " scat2 ", iLabel, s
220 if iLabel == abs(seedLabel):
221 clSeed = np.linalg.inv(clCov)
222 locSeed = np.dot(proL2c, np.dot(clSeed, proL2c.T))
223
224# scatter a little
225 scatNorm = np.random.normal(0., 1., 2)
226 clPar[1:3] = clPar[1:3] + scatErr * scatNorm
227 scatCov[0, 0] += scatErr ** 2
228 scatCov[1, 1] += scatErr ** 2
229# propagate to next measurement layer
230 clPar = np.dot(jacPointToPoint, clPar)
231 clCov = np.dot(jacPointToPoint, np.dot(clCov, jacPointToPoint.T))
232 s += step
233 # propagate scatCov
234 scatCov = np.dot (jacPointToPoint[1:, 1:], np.dot(scatCov, jacPointToPoint[1:, 1:].T))
235 oldL2c = proL2c
236
237# add external seed
238 if locSeed is not None:
239 traj.addExternalSeed(seedLabel, locSeed)
240
241# fit trajectory
242 Chi2, Ndf, Lost, bandCond = traj.fit()
243 print " Record, Chi2, Ndf, Lost", iTry, Chi2, Ndf, Lost, bandCond if bandCond > 0. else 0.
244# dump trajectory
245 #traj.printPoints()
246 #traj.printData()
247# write to MP binary file
248 traj.milleOut(binaryFile)
249# sum up
250 Chi2Sum += Chi2
251 NdfSum += Ndf
252 LostSum += Lost
253# get corrections and covariance matrix at points
254 if (iTry == 0):
255 for i in range(1, nLayer + 1):
256 locPar, locCov = traj.getResults(-i)
257 print " >Point ", i
258 print " locPar ", locPar
259 #print " locCov ", locCov
260 locPar, locCov = traj.getResults(i)
261 print " Point> ", i
262 print " locPar ", locPar
263 #print " locCov ", locCov
264# check residuals
265 for i in range(traj.getNumPoints()):
266 numData, aResiduals, aMeasErr, aResErr, aDownWeight = traj.getMeasResults(i + 1)
267 for j in range(numData):
268 print " measRes " , i, j, aResiduals[j], aMeasErr[j], aResErr[j], aDownWeight[j]
269 #numData, aResiduals, aMeasErr, aResErr, aDownWeight = traj.getScatResults(i + 1)
270 #for j in range(numData):
271 # print " scatRes " , i, j, aResiduals[j], aMeasErr[j], aResErr[j], aDownWeight[j]
272#
273 end = time.clock()
274 print " Time [s] ", end - start
275 print " Chi2Sum/NdfSum ", Chi2Sum / NdfSum
276 print " LostSum/nTry ", LostSum / nTry
277
278
279
282#
283 #binaryFile = open("milleBinaryISN.dat", "rb")
284 #binaryFile = open("/Users/kleinwrt/Documents/workspace/CMSAL/mp3914/milleBinary-151-252458.dat")
285 binaryFile = open("/Users/kleinwrt/Documents/workspace/GBL-thickScat/milleBinaryISN.dat")
286 nRec = 0
287 maxRec = 10 #: maximum number of records to read
288 Chi2Sum = 0.
289 NdfSum = 0
290 LostSum = 0.
291 start = time.clock()
292
293 try:
294 while(nRec < maxRec):
295# create trajectory
296 traj = GblTrajectory(0)
297# read from file
298 traj.milleIn(binaryFile) # get data blocks from file
299 nRec += 1
300# fit trajectory
301 Chi2, Ndf, Lost, bandCond = traj.fit()
302 print " Record, Chi2, Ndf, Lost", nRec, Chi2, Ndf, Lost, bandCond
303# sum up
304 Chi2Sum += Chi2
305 NdfSum += Ndf
306 LostSum += Lost
307
308 except EOFError:
309 pass
310
311 print " records read ", nRec
312 end = time.clock()
313 print " Time [s] ", end - start
314 print " Chi2Sum/NdfSum ", Chi2Sum / NdfSum
315 print " LostSum/nTry ", LostSum / nRec
316
317
318
328def gblSimpleJacobian(ds, cosl, bfac):
329 jac = np.eye(5)
330 jac[1, 0] = -bfac * ds * cosl
331 jac[3, 0] = -0.5 * bfac * ds * ds * cosl
332 jac[3, 1] = ds
333 jac[4, 2] = ds
334 return jac
335
336
337
341class gblMeasSystem(object):
342
343
348 def __init__(self, measDir, curviDir):
349
350 self.__prod = np.dot(curviDir, measDir.T)
351
352
357 return None
358
359
364 meas2crv = np.zeros((5, 5))
365 meas2crv[0, 0] = 1.
366 meas2crv[1:3, 1:3] = self.__prod[1:3, 1:3] * self.__prod[0, 0] # (U,V)*(J,K) * T*I
367 meas2crv[3:5, 3:5] = self.__prod[1:3, 1:3] # (U,V)*(J,K)
368 return meas2crv
369
370
375 return np.linalg.inv(self.__prod[1:3, 1:3]) #((U,V)*(J,K))^-1
376
377
381 def getScatPrecision(self, scatErr):
382 c1 = self.__prod[0, 1] # T*J
383 c2 = self.__prod[0, 2] # T*K
384 fac = (1 - c1 * c1 - c2 * c2) / (scatErr * scatErr)
385 scatP = np.empty((2, 2))
386 scatP[0, 0] = fac * (1 - c1 * c1)
387 scatP[0, 1] = fac * (-c1 * c2)
388 scatP[1, 0] = fac * (-c1 * c2)
389 scatP[1, 1] = fac * (1 - c2 * c2)
390 return scatP
391
392
393
398class gblCurviSystem(object):
399
400
405 def __init__(self, measDir, curviDir):
406
407 self.__c2m = np.linalg.inv(np.dot(curviDir[1:3, 1:3], measDir[1:3, 1:3].T))
408
409
414 return self.__c2m #((U,V)*(J,K))^-1
415
416
421 return np.eye(5) # unit matrix
422
423
428 return self.__c2m #((U,V)*(J,K))^-1
429
430
434 def getScatPrecision(self, scatErr):
435 return np.array([1., 1.]) / (scatErr * scatErr) # diagonal only
436
437
438# create points on initial trajectory, create trajectory from points,
439# fit and write trajectory to MP-II binary file
440# get track parameter corrections and covariance matrix at points
441#example1()
442# read trajectory from MP-II binary file and refit
443example2()
Curvilinear system as local system.
Definition: gbltst.py:398
__c2m
projection curvilinear to measurement system: ((U,V)*(J,K))^-1
Definition: gbltst.py:407
def getTransLocalToCurvi(self)
Transformation of (q/p, slopes, offsets) from local to curvilinear system.
Definition: gbltst.py:420
def __init__(self, measDir, curviDir)
Construct local system.
Definition: gbltst.py:405
def getTransLocalToMeas(self)
Transformation from local to measurement system.
Definition: gbltst.py:413
def getTransCurviToMeas(self)
Transformation from curvilinear to measurement system.
Definition: gbltst.py:427
def getScatPrecision(self, scatErr)
Scattering precision matrix in local system
Definition: gbltst.py:434
Measurement system as local system.
Definition: gbltst.py:341
def __init__(self, measDir, curviDir)
Construct local system.
Definition: gbltst.py:348
def getTransLocalToCurvi(self)
Transformation of (q/p, slopes, offsets) from local to curvilinear system.
Definition: gbltst.py:363
def getTransCurviToMeas(self)
Transformation from curvilinear to measurement system.
Definition: gbltst.py:374
__prod
products (T,U,V) * (I,J,K)
Definition: gbltst.py:350
def getScatPrecision(self, scatErr)
Scattering precision matrix in local system
Definition: gbltst.py:381
def getTransLocalToMeas(self)
Transformation from local to measurement system.
Definition: gbltst.py:356
def example1()
Create points on initial trajectory, create trajectory from points, fit and write trajectory to MP-II...
Definition: gbltst.py:55
def gblSimpleJacobian(ds, cosl, bfac)
Simple jacobian.
Definition: gbltst.py:328
def example2()
Read trajectory from MP-II binary file and refit.
Definition: gbltst.py:281