123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202 |
- from numpy.testing import assert_array_almost_equal, assert_array_equal
- from pytest import raises as assert_raises
- from numpy import array, transpose, dot, conjugate, zeros_like, empty
- from numpy.random import random
- from scipy.linalg import cholesky, cholesky_banded, cho_solve_banded, \
- cho_factor, cho_solve
- from scipy.linalg._testutils import assert_no_overwrite
- class TestCholesky:
- def test_simple(self):
- a = [[8, 2, 3], [2, 9, 3], [3, 3, 6]]
- c = cholesky(a)
- assert_array_almost_equal(dot(transpose(c), c), a)
- c = transpose(c)
- a = dot(c, transpose(c))
- assert_array_almost_equal(cholesky(a, lower=1), c)
- def test_check_finite(self):
- a = [[8, 2, 3], [2, 9, 3], [3, 3, 6]]
- c = cholesky(a, check_finite=False)
- assert_array_almost_equal(dot(transpose(c), c), a)
- c = transpose(c)
- a = dot(c, transpose(c))
- assert_array_almost_equal(cholesky(a, lower=1, check_finite=False), c)
- def test_simple_complex(self):
- m = array([[3+1j, 3+4j, 5], [0, 2+2j, 2+7j], [0, 0, 7+4j]])
- a = dot(transpose(conjugate(m)), m)
- c = cholesky(a)
- a1 = dot(transpose(conjugate(c)), c)
- assert_array_almost_equal(a, a1)
- c = transpose(c)
- a = dot(c, transpose(conjugate(c)))
- assert_array_almost_equal(cholesky(a, lower=1), c)
- def test_random(self):
- n = 20
- for k in range(2):
- m = random([n, n])
- for i in range(n):
- m[i, i] = 20*(.1+m[i, i])
- a = dot(transpose(m), m)
- c = cholesky(a)
- a1 = dot(transpose(c), c)
- assert_array_almost_equal(a, a1)
- c = transpose(c)
- a = dot(c, transpose(c))
- assert_array_almost_equal(cholesky(a, lower=1), c)
- def test_random_complex(self):
- n = 20
- for k in range(2):
- m = random([n, n])+1j*random([n, n])
- for i in range(n):
- m[i, i] = 20*(.1+abs(m[i, i]))
- a = dot(transpose(conjugate(m)), m)
- c = cholesky(a)
- a1 = dot(transpose(conjugate(c)), c)
- assert_array_almost_equal(a, a1)
- c = transpose(c)
- a = dot(c, transpose(conjugate(c)))
- assert_array_almost_equal(cholesky(a, lower=1), c)
- class TestCholeskyBanded:
- """Tests for cholesky_banded() and cho_solve_banded."""
- def test_check_finite(self):
- # Symmetric positive definite banded matrix `a`
- a = array([[4.0, 1.0, 0.0, 0.0],
- [1.0, 4.0, 0.5, 0.0],
- [0.0, 0.5, 4.0, 0.2],
- [0.0, 0.0, 0.2, 4.0]])
- # Banded storage form of `a`.
- ab = array([[-1.0, 1.0, 0.5, 0.2],
- [4.0, 4.0, 4.0, 4.0]])
- c = cholesky_banded(ab, lower=False, check_finite=False)
- ufac = zeros_like(a)
- ufac[list(range(4)), list(range(4))] = c[-1]
- ufac[(0, 1, 2), (1, 2, 3)] = c[0, 1:]
- assert_array_almost_equal(a, dot(ufac.T, ufac))
- b = array([0.0, 0.5, 4.2, 4.2])
- x = cho_solve_banded((c, False), b, check_finite=False)
- assert_array_almost_equal(x, [0.0, 0.0, 1.0, 1.0])
- def test_upper_real(self):
- # Symmetric positive definite banded matrix `a`
- a = array([[4.0, 1.0, 0.0, 0.0],
- [1.0, 4.0, 0.5, 0.0],
- [0.0, 0.5, 4.0, 0.2],
- [0.0, 0.0, 0.2, 4.0]])
- # Banded storage form of `a`.
- ab = array([[-1.0, 1.0, 0.5, 0.2],
- [4.0, 4.0, 4.0, 4.0]])
- c = cholesky_banded(ab, lower=False)
- ufac = zeros_like(a)
- ufac[list(range(4)), list(range(4))] = c[-1]
- ufac[(0, 1, 2), (1, 2, 3)] = c[0, 1:]
- assert_array_almost_equal(a, dot(ufac.T, ufac))
- b = array([0.0, 0.5, 4.2, 4.2])
- x = cho_solve_banded((c, False), b)
- assert_array_almost_equal(x, [0.0, 0.0, 1.0, 1.0])
- def test_upper_complex(self):
- # Hermitian positive definite banded matrix `a`
- a = array([[4.0, 1.0, 0.0, 0.0],
- [1.0, 4.0, 0.5, 0.0],
- [0.0, 0.5, 4.0, -0.2j],
- [0.0, 0.0, 0.2j, 4.0]])
- # Banded storage form of `a`.
- ab = array([[-1.0, 1.0, 0.5, -0.2j],
- [4.0, 4.0, 4.0, 4.0]])
- c = cholesky_banded(ab, lower=False)
- ufac = zeros_like(a)
- ufac[list(range(4)), list(range(4))] = c[-1]
- ufac[(0, 1, 2), (1, 2, 3)] = c[0, 1:]
- assert_array_almost_equal(a, dot(ufac.conj().T, ufac))
- b = array([0.0, 0.5, 4.0-0.2j, 0.2j + 4.0])
- x = cho_solve_banded((c, False), b)
- assert_array_almost_equal(x, [0.0, 0.0, 1.0, 1.0])
- def test_lower_real(self):
- # Symmetric positive definite banded matrix `a`
- a = array([[4.0, 1.0, 0.0, 0.0],
- [1.0, 4.0, 0.5, 0.0],
- [0.0, 0.5, 4.0, 0.2],
- [0.0, 0.0, 0.2, 4.0]])
- # Banded storage form of `a`.
- ab = array([[4.0, 4.0, 4.0, 4.0],
- [1.0, 0.5, 0.2, -1.0]])
- c = cholesky_banded(ab, lower=True)
- lfac = zeros_like(a)
- lfac[list(range(4)), list(range(4))] = c[0]
- lfac[(1, 2, 3), (0, 1, 2)] = c[1, :3]
- assert_array_almost_equal(a, dot(lfac, lfac.T))
- b = array([0.0, 0.5, 4.2, 4.2])
- x = cho_solve_banded((c, True), b)
- assert_array_almost_equal(x, [0.0, 0.0, 1.0, 1.0])
- def test_lower_complex(self):
- # Hermitian positive definite banded matrix `a`
- a = array([[4.0, 1.0, 0.0, 0.0],
- [1.0, 4.0, 0.5, 0.0],
- [0.0, 0.5, 4.0, -0.2j],
- [0.0, 0.0, 0.2j, 4.0]])
- # Banded storage form of `a`.
- ab = array([[4.0, 4.0, 4.0, 4.0],
- [1.0, 0.5, 0.2j, -1.0]])
- c = cholesky_banded(ab, lower=True)
- lfac = zeros_like(a)
- lfac[list(range(4)), list(range(4))] = c[0]
- lfac[(1, 2, 3), (0, 1, 2)] = c[1, :3]
- assert_array_almost_equal(a, dot(lfac, lfac.conj().T))
- b = array([0.0, 0.5j, 3.8j, 3.8])
- x = cho_solve_banded((c, True), b)
- assert_array_almost_equal(x, [0.0, 0.0, 1.0j, 1.0])
- class TestOverwrite:
- def test_cholesky(self):
- assert_no_overwrite(cholesky, [(3, 3)])
- def test_cho_factor(self):
- assert_no_overwrite(cho_factor, [(3, 3)])
- def test_cho_solve(self):
- x = array([[2, -1, 0], [-1, 2, -1], [0, -1, 2]])
- xcho = cho_factor(x)
- assert_no_overwrite(lambda b: cho_solve(xcho, b), [(3,)])
- def test_cholesky_banded(self):
- assert_no_overwrite(cholesky_banded, [(2, 3)])
- def test_cho_solve_banded(self):
- x = array([[0, -1, -1], [2, 2, 2]])
- xcho = cholesky_banded(x)
- assert_no_overwrite(lambda b: cho_solve_banded((xcho, False), b),
- [(3,)])
- class TestEmptyArray:
- def test_cho_factor_empty_square(self):
- a = empty((0, 0))
- b = array([])
- c = array([[]])
- d = []
- e = [[]]
- x, _ = cho_factor(a)
- assert_array_equal(x, a)
- for x in ([b, c, d, e]):
- assert_raises(ValueError, cho_factor, x)
|