1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429
|
#!/usr/bin/env python
#
# Created by: Pearu Peterson, March 2002
#
""" Test functions for linalg.basic module
"""
"""
Bugs:
1) solve.check_random_sym_complex fails if a is complex
and transpose(a) = conjugate(a) (a is Hermitian).
"""
__usage__ = """
Build linalg:
python setup_linalg.py build
Run tests if scipy is installed:
python -c 'import scipy;scipy.linalg.test()'
Run tests if linalg is not installed:
python tests/test_basic.py
"""
from numpy import arange, add, array, dot, zeros, identity, conjugate, transpose
import numpy.linalg as linalg
from numpy.testing import *
from scipy.linalg import solve,inv,det,lstsq, toeplitz, hankel, tri, triu, \
tril, pinv, pinv2, solve_banded
def random(size):
return rand(*size)
def get_mat(n):
data = arange(n)
data = add.outer(data,data)
return data
class TestSolveBanded(TestCase):
def test_simple(self):
a = [[1,20,0,0],[-30,4,6,0],[2,1,20,2],[0,-1,7,14]]
ab = [[0,20,6,2],
[1,4,20,14],
[-30,1,7,0],
[2,-1,0,0]]
l,u = 2,1
for b in ([[1,0,0,0],[0,0,0,1],[0,1,0,0],[0,1,0,0]],
[[2,1],[-30,4],[2,3],[1,3]]):
x = solve_banded((l,u),ab,b)
assert_array_almost_equal(dot(a,x),b)
class TestSolve(TestCase):
def test_20Feb04_bug(self):
a = [[1,1],[1.0,0]] # ok
x0 = solve(a,[1,0j])
assert_array_almost_equal(dot(a,x0),[1,0])
a = [[1,1],[1.2,0]] # gives failure with clapack.zgesv(..,rowmajor=0)
b = [1,0j]
x0 = solve(a,b)
assert_array_almost_equal(dot(a,x0),[1,0])
def test_simple(self):
a = [[1,20],[-30,4]]
for b in ([[1,0],[0,1]],[1,0],
[[2,1],[-30,4]]):
x = solve(a,b)
assert_array_almost_equal(dot(a,x),b)
def test_simple_sym(self):
a = [[2,3],[3,5]]
for lower in [0,1]:
for b in ([[1,0],[0,1]],[1,0]):
x = solve(a,b,sym_pos=1,lower=lower)
assert_array_almost_equal(dot(a,x),b)
def test_simple_sym_complex(self):
a = [[5,2],[2,4]]
for b in [[1j,0],
[[1j,1j],
[0,2]],
]:
x = solve(a,b,sym_pos=1)
assert_array_almost_equal(dot(a,x),b)
def test_simple_complex(self):
a = array([[5,2],[2j,4]],'D')
for b in [[1j,0],
[[1j,1j],
[0,2]],
[1,0j],
array([1,0],'D'),
]:
x = solve(a,b)
assert_array_almost_equal(dot(a,x),b)
def test_nils_20Feb04(self):
n = 2
A = random([n,n])+random([n,n])*1j
X = zeros((n,n),'D')
Ainv = inv(A)
R = identity(n)+identity(n)*0j
for i in arange(0,n):
r = R[:,i]
X[:,i] = solve(A,r)
assert_array_almost_equal(X,Ainv)
def test_random(self):
n = 20
a = random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
for i in range(4):
b = random([n,3])
x = solve(a,b)
assert_array_almost_equal(dot(a,x),b)
def test_random_complex(self):
n = 20
a = random([n,n]) + 1j * random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
for i in range(2):
b = random([n,3])
x = solve(a,b)
assert_array_almost_equal(dot(a,x),b)
def test_random_sym(self):
n = 20
a = random([n,n])
for i in range(n):
a[i,i] = abs(20*(.1+a[i,i]))
for j in range(i):
a[i,j] = a[j,i]
for i in range(4):
b = random([n])
x = solve(a,b,sym_pos=1)
assert_array_almost_equal(dot(a,x),b)
def test_random_sym_complex(self):
n = 20
a = random([n,n])
#a = a + 1j*random([n,n]) # XXX: with this the accuracy will be very low
for i in range(n):
a[i,i] = abs(20*(.1+a[i,i]))
for j in range(i):
a[i,j] = conjugate(a[j,i])
b = random([n])+2j*random([n])
for i in range(2):
x = solve(a,b,sym_pos=1)
assert_array_almost_equal(dot(a,x),b)
class TestInv(TestCase):
def test_simple(self):
a = [[1,2],[3,4]]
a_inv = inv(a)
assert_array_almost_equal(dot(a,a_inv),
[[1,0],[0,1]])
a = [[1,2,3],[4,5,6],[7,8,10]]
a_inv = inv(a)
assert_array_almost_equal(dot(a,a_inv),
[[1,0,0],[0,1,0],[0,0,1]])
def test_random(self):
n = 20
for i in range(4):
a = random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
a_inv = inv(a)
assert_array_almost_equal(dot(a,a_inv),
identity(n))
def test_simple_complex(self):
a = [[1,2],[3,4j]]
a_inv = inv(a)
assert_array_almost_equal(dot(a,a_inv),
[[1,0],[0,1]])
def test_random_complex(self):
n = 20
for i in range(4):
a = random([n,n])+2j*random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
a_inv = inv(a)
assert_array_almost_equal(dot(a,a_inv),
identity(n))
class TestDet(TestCase):
def test_simple(self):
a = [[1,2],[3,4]]
a_det = det(a)
assert_almost_equal(a_det,-2.0)
def test_simple_complex(self):
a = [[1,2],[3,4j]]
a_det = det(a)
assert_almost_equal(a_det,-6+4j)
def test_random(self):
basic_det = linalg.det
n = 20
for i in range(4):
a = random([n,n])
d1 = det(a)
d2 = basic_det(a)
assert_almost_equal(d1,d2)
def test_random_complex(self):
basic_det = linalg.det
n = 20
for i in range(4):
a = random([n,n]) + 2j*random([n,n])
d1 = det(a)
d2 = basic_det(a)
assert_almost_equal(d1,d2)
def direct_lstsq(a,b,cmplx=0):
at = transpose(a)
if cmplx:
at = conjugate(at)
a1 = dot(at, a)
b1 = dot(at, b)
return solve(a1, b1)
class TestLstsq(TestCase):
def test_random_overdet_large(self):
#bug report: Nils Wagner
n = 200
a = random([n,2])
for i in range(2): a[i,i] = 20*(.1+a[i,i])
b = random([n,3])
x = lstsq(a,b)[0]
assert_array_almost_equal(x,direct_lstsq(a,b))
def test_simple_exact(self):
a = [[1,20],[-30,4]]
for b in ([[1,0],[0,1]],[1,0],
[[2,1],[-30,4]]):
x = lstsq(a,b)[0]
assert_array_almost_equal(dot(a,x),b)
def test_simple_overdet(self):
a = [[1,2],[4,5],[3,4]]
b = [1,2,3]
x,res,r,s = lstsq(a,b)
#XXX: check defintion of res
assert_array_almost_equal(x,direct_lstsq(a,b))
def test_simple_underdet(self):
a = [[1,2,3],[4,5,6]]
b = [1,2]
x,res,r,s = lstsq(a,b)
#XXX: need independent check
assert_array_almost_equal(x,[[-0.05555556],
[0.11111111],[0.27777778]])
def test_random_exact(self):
n = 20
a = random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
for i in range(4):
b = random([n,3])
x = lstsq(a,b)[0]
assert_array_almost_equal(dot(a,x),b)
def test_random_complex_exact(self):
n = 20
a = random([n,n]) + 1j * random([n,n])
for i in range(n): a[i,i] = 20*(.1+a[i,i])
for i in range(2):
b = random([n,3])
x = lstsq(a,b)[0]
assert_array_almost_equal(dot(a,x),b)
def test_random_overdet(self):
n = 20
m = 15
a = random([n,m])
for i in range(m): a[i,i] = 20*(.1+a[i,i])
for i in range(4):
b = random([n,3])
x,res,r,s = lstsq(a,b)
assert r==m,'unexpected efficient rank'
#XXX: check definition of res
assert_array_almost_equal(x,direct_lstsq(a,b))
def test_random_complex_overdet(self):
n = 20
m = 15
a = random([n,m]) + 1j * random([n,m])
for i in range(m):
a[i,i] = 20*(.1+a[i,i])
for i in range(2):
b = random([n,3])
x,res,r,s = lstsq(a,b)
assert r==m,'unexpected efficient rank'
#XXX: check definition of res
assert_array_almost_equal(x,direct_lstsq(a,b,1))
class TestTri(TestCase):
def test_basic(self):
assert_equal(tri(4),array([[1,0,0,0],
[1,1,0,0],
[1,1,1,0],
[1,1,1,1]]))
assert_equal(tri(4,dtype='f'),array([[1,0,0,0],
[1,1,0,0],
[1,1,1,0],
[1,1,1,1]],'f'))
def test_diag(self):
assert_equal(tri(4,k=1),array([[1,1,0,0],
[1,1,1,0],
[1,1,1,1],
[1,1,1,1]]))
assert_equal(tri(4,k=-1),array([[0,0,0,0],
[1,0,0,0],
[1,1,0,0],
[1,1,1,0]]))
def test_2d(self):
assert_equal(tri(4,3),array([[1,0,0],
[1,1,0],
[1,1,1],
[1,1,1]]))
assert_equal(tri(3,4),array([[1,0,0,0],
[1,1,0,0],
[1,1,1,0]]))
def test_diag2d(self):
assert_equal(tri(3,4,k=2),array([[1,1,1,0],
[1,1,1,1],
[1,1,1,1]]))
assert_equal(tri(4,3,k=-2),array([[0,0,0],
[0,0,0],
[1,0,0],
[1,1,0]]))
class TestTril(TestCase):
def test_basic(self):
a = (100*get_mat(5)).astype('l')
b = a.copy()
for k in range(5):
for l in range(k+1,5):
b[k,l] = 0
assert_equal(tril(a),b)
def test_diag(self):
a = (100*get_mat(5)).astype('f')
b = a.copy()
for k in range(5):
for l in range(k+3,5):
b[k,l] = 0
assert_equal(tril(a,k=2),b)
b = a.copy()
for k in range(5):
for l in range(max((k-1,0)),5):
b[k,l] = 0
assert_equal(tril(a,k=-2),b)
class TestTriu(TestCase):
def test_basic(self):
a = (100*get_mat(5)).astype('l')
b = a.copy()
for k in range(5):
for l in range(k+1,5):
b[l,k] = 0
assert_equal(triu(a),b)
def test_diag(self):
a = (100*get_mat(5)).astype('f')
b = a.copy()
for k in range(5):
for l in range(max((k-1,0)),5):
b[l,k] = 0
assert_equal(triu(a,k=2),b)
b = a.copy()
for k in range(5):
for l in range(k+3,5):
b[l,k] = 0
assert_equal(triu(a,k=-2),b)
class TestToeplitz(TestCase):
def test_basic(self):
y = toeplitz([1,2,3])
assert_array_equal(y,[[1,2,3],[2,1,2],[3,2,1]])
y = toeplitz([1,2,3],[1,4,5])
assert_array_equal(y,[[1,4,5],[2,1,4],[3,2,1]])
class TestHankel(TestCase):
def test_basic(self):
y = hankel([1,2,3])
assert_array_equal(y,[[1,2,3],[2,3,0],[3,0,0]])
y = hankel([1,2,3],[3,4,5])
assert_array_equal(y,[[1,2,3],[2,3,4],[3,4,5]])
class TestPinv(TestCase):
def test_simple(self):
a=array([[1,2,3],[4,5,6.],[7,8,10]])
a_pinv = pinv(a)
assert_array_almost_equal(dot(a,a_pinv),[[1,0,0],[0,1,0],[0,0,1]])
a_pinv = pinv2(a)
assert_array_almost_equal(dot(a,a_pinv),[[1,0,0],[0,1,0],[0,0,1]])
def test_simple_0det(self):
a=array([[1,2,3],[4,5,6.],[7,8,9]])
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
def test_simple_cols(self):
a=array([[1,2,3],[4,5,6.]])
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
def test_simple_rows(self):
a=array([[1,2],[3,4],[5,6]])
a_pinv = pinv(a)
a_pinv2 = pinv2(a)
assert_array_almost_equal(a_pinv,a_pinv2)
if __name__ == "__main__":
run_module_suite()
|