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 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526
|
#!/usr/bin/env python
#
# Created by: Pearu Peterson, March 2002
#
""" Test functions for scipy.linalg.matfuncs module
"""
from __future__ import division, print_function, absolute_import
import math
import warnings
import numpy as np
from numpy import array, eye, exp, random
from numpy.linalg import matrix_power
from numpy.testing import (TestCase, run_module_suite,
assert_allclose, assert_, assert_array_almost_equal, assert_equal,
assert_array_almost_equal_nulp)
from scipy.sparse import csc_matrix, SparseEfficiencyWarning
from scipy.sparse.construct import eye as speye
from scipy.sparse.linalg.matfuncs import (expm, _expm,
ProductOperator, MatrixPowerOperator,
_onenorm_matrix_power_nnm)
from scipy.linalg import logm
from scipy.special import factorial
import scipy.sparse
import scipy.sparse.linalg
def _burkardt_13_power(n, p):
"""
A helper function for testing matrix functions.
Parameters
----------
n : integer greater than 1
Order of the square matrix to be returned.
p : non-negative integer
Power of the matrix.
Returns
-------
out : ndarray representing a square matrix
A Forsythe matrix of order n, raised to the power p.
"""
# Input validation.
if n != int(n) or n < 2:
raise ValueError('n must be an integer greater than 1')
n = int(n)
if p != int(p) or p < 0:
raise ValueError('p must be a non-negative integer')
p = int(p)
# Construct the matrix explicitly.
a, b = divmod(p, n)
large = np.power(10.0, -n*a)
small = large * np.power(10.0, -n)
return np.diag([large]*(n-b), b) + np.diag([small]*b, b-n)
def test_onenorm_matrix_power_nnm():
np.random.seed(1234)
for n in range(1, 5):
for p in range(5):
M = np.random.random((n, n))
Mp = np.linalg.matrix_power(M, p)
observed = _onenorm_matrix_power_nnm(M, p)
expected = np.linalg.norm(Mp, 1)
assert_allclose(observed, expected)
class TestExpM(TestCase):
def test_zero_ndarray(self):
a = array([[0.,0],[0,0]])
assert_array_almost_equal(expm(a),[[1,0],[0,1]])
def test_zero_sparse(self):
a = csc_matrix([[0.,0],[0,0]])
assert_array_almost_equal(expm(a).toarray(),[[1,0],[0,1]])
def test_zero_matrix(self):
a = np.matrix([[0.,0],[0,0]])
assert_array_almost_equal(expm(a),[[1,0],[0,1]])
def test_misc_types(self):
A = expm(np.array([[1]]))
yield assert_allclose, expm(((1,),)), A
yield assert_allclose, expm([[1]]), A
yield assert_allclose, expm(np.matrix([[1]])), A
yield assert_allclose, expm(np.array([[1]])), A
yield assert_allclose, expm(csc_matrix([[1]])), A
B = expm(np.array([[1j]]))
yield assert_allclose, expm(((1j,),)), B
yield assert_allclose, expm([[1j]]), B
yield assert_allclose, expm(np.matrix([[1j]])), B
yield assert_allclose, expm(csc_matrix([[1j]])), B
def test_bidiagonal_sparse(self):
A = csc_matrix([
[1, 3, 0],
[0, 1, 5],
[0, 0, 2]], dtype=float)
e1 = math.exp(1)
e2 = math.exp(2)
expected = np.array([
[e1, 3*e1, 15*(e2 - 2*e1)],
[0, e1, 5*(e2 - e1)],
[0, 0, e2]], dtype=float)
observed = expm(A).toarray()
assert_array_almost_equal(observed, expected)
def test_padecases_dtype_float(self):
for dtype in [np.float32, np.float64]:
for scale in [1e-2, 1e-1, 5e-1, 1, 10]:
A = scale * eye(3, dtype=dtype)
observed = expm(A)
expected = exp(scale) * eye(3, dtype=dtype)
assert_array_almost_equal_nulp(observed, expected, nulp=100)
def test_padecases_dtype_complex(self):
for dtype in [np.complex64, np.complex128]:
for scale in [1e-2, 1e-1, 5e-1, 1, 10]:
A = scale * eye(3, dtype=dtype)
observed = expm(A)
expected = exp(scale) * eye(3, dtype=dtype)
assert_array_almost_equal_nulp(observed, expected, nulp=100)
def test_padecases_dtype_sparse_float(self):
# float32 and complex64 lead to errors in spsolve/UMFpack
dtype = np.float64
with warnings.catch_warnings():
warnings.simplefilter("ignore", category=SparseEfficiencyWarning)
for scale in [1e-2, 1e-1, 5e-1, 1, 10]:
a = scale * speye(3, 3, dtype=dtype, format='csc')
e = exp(scale) * eye(3, dtype=dtype)
exact_onenorm = _expm(a, use_exact_onenorm=True).toarray()
inexact_onenorm = _expm(a, use_exact_onenorm=False).toarray()
assert_array_almost_equal_nulp(exact_onenorm, e, nulp=100)
assert_array_almost_equal_nulp(inexact_onenorm, e, nulp=100)
def test_padecases_dtype_sparse_complex(self):
# float32 and complex64 lead to errors in spsolve/UMFpack
dtype = np.complex128
with warnings.catch_warnings():
warnings.simplefilter("ignore", category=SparseEfficiencyWarning)
for scale in [1e-2, 1e-1, 5e-1, 1, 10]:
a = scale * speye(3, 3, dtype=dtype, format='csc')
e = exp(scale) * eye(3, dtype=dtype)
assert_array_almost_equal_nulp(expm(a).toarray(), e, nulp=100)
def test_logm_consistency(self):
random.seed(1234)
for dtype in [np.float64, np.complex128]:
for n in range(1, 10):
for scale in [1e-4, 1e-3, 1e-2, 1e-1, 1, 1e1, 1e2]:
# make logm(A) be of a given scale
A = (eye(n) + random.rand(n, n) * scale).astype(dtype)
if np.iscomplexobj(A):
A = A + 1j * random.rand(n, n) * scale
assert_array_almost_equal(expm(logm(A)), A)
def test_integer_matrix(self):
Q = np.array([
[-3, 1, 1, 1],
[1, -3, 1, 1],
[1, 1, -3, 1],
[1, 1, 1, -3]])
assert_allclose(expm(Q), expm(1.0 * Q))
def test_triangularity_perturbation(self):
# Experiment (1) of
# Awad H. Al-Mohy and Nicholas J. Higham (2012)
# Improved Inverse Scaling and Squaring Algorithms
# for the Matrix Logarithm.
A = np.array([
[3.2346e-1, 3e4, 3e4, 3e4],
[0, 3.0089e-1, 3e4, 3e4],
[0, 0, 3.221e-1, 3e4],
[0, 0, 0, 3.0744e-1]],
dtype=float)
A_logm = np.array([
[-1.12867982029050462e+00, 9.61418377142025565e+04,
-4.52485573953179264e+09, 2.92496941103871812e+14],
[0.00000000000000000e+00, -1.20101052953082288e+00,
9.63469687211303099e+04, -4.68104828911105442e+09],
[0.00000000000000000e+00, 0.00000000000000000e+00,
-1.13289322264498393e+00, 9.53249183094775653e+04],
[0.00000000000000000e+00, 0.00000000000000000e+00,
0.00000000000000000e+00, -1.17947533272554850e+00]],
dtype=float)
assert_allclose(expm(A_logm), A, rtol=1e-4)
# Perturb the upper triangular matrix by tiny amounts,
# so that it becomes technically not upper triangular.
random.seed(1234)
tiny = 1e-17
A_logm_perturbed = A_logm.copy()
A_logm_perturbed[1, 0] = tiny
A_expm_logm_perturbed = expm(A_logm_perturbed)
rtol = 1e-4
atol = 100 * tiny
assert_(not np.allclose(A_expm_logm_perturbed, A, rtol=rtol, atol=atol))
def test_burkardt_1(self):
# This matrix is diagonal.
# The calculation of the matrix exponential is simple.
#
# This is the first of a series of matrix exponential tests
# collected by John Burkardt from the following sources.
#
# Alan Laub,
# Review of "Linear System Theory" by Joao Hespanha,
# SIAM Review,
# Volume 52, Number 4, December 2010, pages 779--781.
#
# Cleve Moler and Charles Van Loan,
# Nineteen Dubious Ways to Compute the Exponential of a Matrix,
# Twenty-Five Years Later,
# SIAM Review,
# Volume 45, Number 1, March 2003, pages 3--49.
#
# Cleve Moler,
# Cleve's Corner: A Balancing Act for the Matrix Exponential,
# 23 July 2012.
#
# Robert Ward,
# Numerical computation of the matrix exponential
# with accuracy estimate,
# SIAM Journal on Numerical Analysis,
# Volume 14, Number 4, September 1977, pages 600--610.
exp1 = np.exp(1)
exp2 = np.exp(2)
A = np.array([
[1, 0],
[0, 2],
], dtype=float)
desired = np.array([
[exp1, 0],
[0, exp2],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_2(self):
# This matrix is symmetric.
# The calculation of the matrix exponential is straightforward.
A = np.array([
[1, 3],
[3, 2],
], dtype=float)
desired = np.array([
[39.322809708033859, 46.166301438885753],
[46.166301438885768, 54.711576854329110],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_3(self):
# This example is due to Laub.
# This matrix is ill-suited for the Taylor series approach.
# As powers of A are computed, the entries blow up too quickly.
exp1 = np.exp(1)
exp39 = np.exp(39)
A = np.array([
[0, 1],
[-39, -40],
], dtype=float)
desired = np.array([
[
39/(38*exp1) - 1/(38*exp39),
-np.expm1(-38) / (38*exp1)],
[
39*np.expm1(-38) / (38*exp1),
-1/(38*exp1) + 39/(38*exp39)],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_4(self):
# This example is due to Moler and Van Loan.
# The example will cause problems for the series summation approach,
# as well as for diagonal Pade approximations.
A = np.array([
[-49, 24],
[-64, 31],
], dtype=float)
U = np.array([[3, 1], [4, 2]], dtype=float)
V = np.array([[1, -1/2], [-2, 3/2]], dtype=float)
w = np.array([-17, -1], dtype=float)
desired = np.dot(U * np.exp(w), V)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_5(self):
# This example is due to Moler and Van Loan.
# This matrix is strictly upper triangular
# All powers of A are zero beyond some (low) limit.
# This example will cause problems for Pade approximations.
A = np.array([
[0, 6, 0, 0],
[0, 0, 6, 0],
[0, 0, 0, 6],
[0, 0, 0, 0],
], dtype=float)
desired = np.array([
[1, 6, 18, 36],
[0, 1, 6, 18],
[0, 0, 1, 6],
[0, 0, 0, 1],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_6(self):
# This example is due to Moler and Van Loan.
# This matrix does not have a complete set of eigenvectors.
# That means the eigenvector approach will fail.
exp1 = np.exp(1)
A = np.array([
[1, 1],
[0, 1],
], dtype=float)
desired = np.array([
[exp1, exp1],
[0, exp1],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_7(self):
# This example is due to Moler and Van Loan.
# This matrix is very close to example 5.
# Mathematically, it has a complete set of eigenvectors.
# Numerically, however, the calculation will be suspect.
exp1 = np.exp(1)
eps = np.spacing(1)
A = np.array([
[1 + eps, 1],
[0, 1 - eps],
], dtype=float)
desired = np.array([
[exp1, exp1],
[0, exp1],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_8(self):
# This matrix was an example in Wikipedia.
exp4 = np.exp(4)
exp16 = np.exp(16)
A = np.array([
[21, 17, 6],
[-5, -1, -6],
[4, 4, 16],
], dtype=float)
desired = np.array([
[13*exp16 - exp4, 13*exp16 - 5*exp4, 2*exp16 - 2*exp4],
[-9*exp16 + exp4, -9*exp16 + 5*exp4, -2*exp16 + 2*exp4],
[16*exp16, 16*exp16, 4*exp16],
], dtype=float) * 0.25
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_9(self):
# This matrix is due to the NAG Library.
# It is an example for function F01ECF.
A = np.array([
[1, 2, 2, 2],
[3, 1, 1, 2],
[3, 2, 1, 2],
[3, 3, 3, 1],
], dtype=float)
desired = np.array([
[740.7038, 610.8500, 542.2743, 549.1753],
[731.2510, 603.5524, 535.0884, 542.2743],
[823.7630, 679.4257, 603.5524, 610.8500],
[998.4355, 823.7630, 731.2510, 740.7038],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_10(self):
# This is Ward's example #1.
# It is defective and nonderogatory.
A = np.array([
[4, 2, 0],
[1, 4, 1],
[1, 1, 4],
], dtype=float)
assert_allclose(sorted(scipy.linalg.eigvals(A)), (3, 3, 6))
desired = np.array([
[147.8666224463699, 183.7651386463682, 71.79703239999647],
[127.7810855231823, 183.7651386463682, 91.88256932318415],
[127.7810855231824, 163.6796017231806, 111.9681062463718],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_11(self):
# This is Ward's example #2.
# It is a symmetric matrix.
A = np.array([
[29.87942128909879, 0.7815750847907159, -2.289519314033932],
[0.7815750847907159, 25.72656945571064, 8.680737820540137],
[-2.289519314033932, 8.680737820540137, 34.39400925519054],
], dtype=float)
assert_allclose(scipy.linalg.eigvalsh(A), (20, 30, 40))
desired = np.array([
[
5.496313853692378E+15,
-1.823188097200898E+16,
-3.047577080858001E+16],
[
-1.823188097200899E+16,
6.060522870222108E+16,
1.012918429302482E+17],
[
-3.047577080858001E+16,
1.012918429302482E+17,
1.692944112408493E+17],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_12(self):
# This is Ward's example #3.
# Ward's algorithm has difficulty estimating the accuracy
# of its results.
A = np.array([
[-131, 19, 18],
[-390, 56, 54],
[-387, 57, 52],
], dtype=float)
assert_allclose(sorted(scipy.linalg.eigvals(A)), (-20, -2, -1))
desired = np.array([
[-1.509644158793135, 0.3678794391096522, 0.1353352811751005],
[-5.632570799891469, 1.471517758499875, 0.4060058435250609],
[-4.934938326088363, 1.103638317328798, 0.5413411267617766],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
def test_burkardt_13(self):
# This is Ward's example #4.
# This is a version of the Forsythe matrix.
# The eigenvector problem is badly conditioned.
# Ward's algorithm has difficulty esimating the accuracy
# of its results for this problem.
#
# Check the construction of one instance of this family of matrices.
A4_actual = _burkardt_13_power(4, 1)
A4_desired = [[0, 1, 0, 0],
[0, 0, 1, 0],
[0, 0, 0, 1],
[1e-4, 0, 0, 0]]
assert_allclose(A4_actual, A4_desired)
# Check the expm for a few instances.
for n in (2, 3, 4, 10):
# Approximate expm using Taylor series.
# This works well for this matrix family
# because each matrix in the summation,
# even before dividing by the factorial,
# is entrywise positive with max entry 10**(-floor(p/n)*n).
k = max(1, int(np.ceil(16/n)))
desired = np.zeros((n, n), dtype=float)
for p in range(n*k):
Ap = _burkardt_13_power(n, p)
assert_equal(np.min(Ap), 0)
assert_allclose(np.max(Ap), np.power(10, -np.floor(p/n)*n))
desired += Ap / factorial(p)
actual = expm(_burkardt_13_power(n, 1))
assert_allclose(actual, desired)
def test_burkardt_14(self):
# This is Moler's example.
# This badly scaled matrix caused problems for MATLAB's expm().
A = np.array([
[0, 1e-8, 0],
[-(2e10 + 4e8/6.), -3, 2e10],
[200./3., 0, -200./3.],
], dtype=float)
desired = np.array([
[0.446849468283175, 1.54044157383952e-09, 0.462811453558774],
[-5743067.77947947, -0.0152830038686819, -4526542.71278401],
[0.447722977849494, 1.54270484519591e-09, 0.463480648837651],
], dtype=float)
actual = expm(A)
assert_allclose(actual, desired)
class TestOperators(TestCase):
def test_product_operator(self):
random.seed(1234)
n = 5
k = 2
nsamples = 10
for i in range(nsamples):
A = np.random.randn(n, n)
B = np.random.randn(n, n)
C = np.random.randn(n, n)
D = np.random.randn(n, k)
op = ProductOperator(A, B, C)
assert_allclose(op.matmat(D), A.dot(B).dot(C).dot(D))
assert_allclose(op.T.matmat(D), (A.dot(B).dot(C)).T.dot(D))
def test_matrix_power_operator(self):
random.seed(1234)
n = 5
k = 2
p = 3
nsamples = 10
for i in range(nsamples):
A = np.random.randn(n, n)
B = np.random.randn(n, k)
op = MatrixPowerOperator(A, p)
assert_allclose(op.matmat(B), matrix_power(A, p).dot(B))
assert_allclose(op.T.matmat(B), matrix_power(A, p).T.dot(B))
if __name__ == "__main__":
run_module_suite()
|