test_matfuncs.py 32 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836
  1. #
  2. # Created by: Pearu Peterson, March 2002
  3. #
  4. """ Test functions for linalg.matfuncs module
  5. """
  6. from __future__ import division, print_function, absolute_import
  7. import random
  8. import functools
  9. import numpy as np
  10. from numpy import array, matrix, identity, dot, sqrt, double
  11. from numpy.testing import (
  12. assert_array_equal, assert_array_less, assert_equal,
  13. assert_array_almost_equal, assert_array_almost_equal_nulp,
  14. assert_allclose, assert_)
  15. import pytest
  16. from scipy._lib._numpy_compat import _assert_warns, suppress_warnings
  17. import scipy.linalg
  18. from scipy.linalg import (funm, signm, logm, sqrtm, fractional_matrix_power,
  19. expm, expm_frechet, expm_cond, norm)
  20. from scipy.linalg import _matfuncs_inv_ssq
  21. import scipy.linalg._expm_frechet
  22. from scipy.optimize import minimize
  23. def _get_al_mohy_higham_2012_experiment_1():
  24. """
  25. Return the test matrix from Experiment (1) of [1]_.
  26. References
  27. ----------
  28. .. [1] Awad H. Al-Mohy and Nicholas J. Higham (2012)
  29. "Improved Inverse Scaling and Squaring Algorithms
  30. for the Matrix Logarithm."
  31. SIAM Journal on Scientific Computing, 34 (4). C152-C169.
  32. ISSN 1095-7197
  33. """
  34. A = np.array([
  35. [3.2346e-1, 3e4, 3e4, 3e4],
  36. [0, 3.0089e-1, 3e4, 3e4],
  37. [0, 0, 3.2210e-1, 3e4],
  38. [0, 0, 0, 3.0744e-1]], dtype=float)
  39. return A
  40. class TestSignM(object):
  41. def test_nils(self):
  42. a = array([[29.2, -24.2, 69.5, 49.8, 7.],
  43. [-9.2, 5.2, -18., -16.8, -2.],
  44. [-10., 6., -20., -18., -2.],
  45. [-9.6, 9.6, -25.5, -15.4, -2.],
  46. [9.8, -4.8, 18., 18.2, 2.]])
  47. cr = array([[11.94933333,-2.24533333,15.31733333,21.65333333,-2.24533333],
  48. [-3.84266667,0.49866667,-4.59066667,-7.18666667,0.49866667],
  49. [-4.08,0.56,-4.92,-7.6,0.56],
  50. [-4.03466667,1.04266667,-5.59866667,-7.02666667,1.04266667],
  51. [4.15733333,-0.50133333,4.90933333,7.81333333,-0.50133333]])
  52. r = signm(a)
  53. assert_array_almost_equal(r,cr)
  54. def test_defective1(self):
  55. a = array([[0.0,1,0,0],[1,0,1,0],[0,0,0,1],[0,0,1,0]])
  56. r = signm(a, disp=False)
  57. #XXX: what would be the correct result?
  58. def test_defective2(self):
  59. a = array((
  60. [29.2,-24.2,69.5,49.8,7.0],
  61. [-9.2,5.2,-18.0,-16.8,-2.0],
  62. [-10.0,6.0,-20.0,-18.0,-2.0],
  63. [-9.6,9.6,-25.5,-15.4,-2.0],
  64. [9.8,-4.8,18.0,18.2,2.0]))
  65. r = signm(a, disp=False)
  66. #XXX: what would be the correct result?
  67. def test_defective3(self):
  68. a = array([[-2., 25., 0., 0., 0., 0., 0.],
  69. [0., -3., 10., 3., 3., 3., 0.],
  70. [0., 0., 2., 15., 3., 3., 0.],
  71. [0., 0., 0., 0., 15., 3., 0.],
  72. [0., 0., 0., 0., 3., 10., 0.],
  73. [0., 0., 0., 0., 0., -2., 25.],
  74. [0., 0., 0., 0., 0., 0., -3.]])
  75. r = signm(a, disp=False)
  76. #XXX: what would be the correct result?
  77. class TestLogM(object):
  78. def test_nils(self):
  79. a = array([[-2., 25., 0., 0., 0., 0., 0.],
  80. [0., -3., 10., 3., 3., 3., 0.],
  81. [0., 0., 2., 15., 3., 3., 0.],
  82. [0., 0., 0., 0., 15., 3., 0.],
  83. [0., 0., 0., 0., 3., 10., 0.],
  84. [0., 0., 0., 0., 0., -2., 25.],
  85. [0., 0., 0., 0., 0., 0., -3.]])
  86. m = (identity(7)*3.1+0j)-a
  87. logm(m, disp=False)
  88. #XXX: what would be the correct result?
  89. def test_al_mohy_higham_2012_experiment_1_logm(self):
  90. # The logm completes the round trip successfully.
  91. # Note that the expm leg of the round trip is badly conditioned.
  92. A = _get_al_mohy_higham_2012_experiment_1()
  93. A_logm, info = logm(A, disp=False)
  94. A_round_trip = expm(A_logm)
  95. assert_allclose(A_round_trip, A, rtol=1e-5, atol=1e-14)
  96. def test_al_mohy_higham_2012_experiment_1_funm_log(self):
  97. # The raw funm with np.log does not complete the round trip.
  98. # Note that the expm leg of the round trip is badly conditioned.
  99. A = _get_al_mohy_higham_2012_experiment_1()
  100. A_funm_log, info = funm(A, np.log, disp=False)
  101. A_round_trip = expm(A_funm_log)
  102. assert_(not np.allclose(A_round_trip, A, rtol=1e-5, atol=1e-14))
  103. def test_round_trip_random_float(self):
  104. np.random.seed(1234)
  105. for n in range(1, 6):
  106. M_unscaled = np.random.randn(n, n)
  107. for scale in np.logspace(-4, 4, 9):
  108. M = M_unscaled * scale
  109. # Eigenvalues are related to the branch cut.
  110. W = np.linalg.eigvals(M)
  111. err_msg = 'M:{0} eivals:{1}'.format(M, W)
  112. # Check sqrtm round trip because it is used within logm.
  113. M_sqrtm, info = sqrtm(M, disp=False)
  114. M_sqrtm_round_trip = M_sqrtm.dot(M_sqrtm)
  115. assert_allclose(M_sqrtm_round_trip, M)
  116. # Check logm round trip.
  117. M_logm, info = logm(M, disp=False)
  118. M_logm_round_trip = expm(M_logm)
  119. assert_allclose(M_logm_round_trip, M, err_msg=err_msg)
  120. def test_round_trip_random_complex(self):
  121. np.random.seed(1234)
  122. for n in range(1, 6):
  123. M_unscaled = np.random.randn(n, n) + 1j * np.random.randn(n, n)
  124. for scale in np.logspace(-4, 4, 9):
  125. M = M_unscaled * scale
  126. M_logm, info = logm(M, disp=False)
  127. M_round_trip = expm(M_logm)
  128. assert_allclose(M_round_trip, M)
  129. def test_logm_type_preservation_and_conversion(self):
  130. # The logm matrix function should preserve the type of a matrix
  131. # whose eigenvalues are positive with zero imaginary part.
  132. # Test this preservation for variously structured matrices.
  133. complex_dtype_chars = ('F', 'D', 'G')
  134. for matrix_as_list in (
  135. [[1, 0], [0, 1]],
  136. [[1, 0], [1, 1]],
  137. [[2, 1], [1, 1]],
  138. [[2, 3], [1, 2]]):
  139. # check that the spectrum has the expected properties
  140. W = scipy.linalg.eigvals(matrix_as_list)
  141. assert_(not any(w.imag or w.real < 0 for w in W))
  142. # check float type preservation
  143. A = np.array(matrix_as_list, dtype=float)
  144. A_logm, info = logm(A, disp=False)
  145. assert_(A_logm.dtype.char not in complex_dtype_chars)
  146. # check complex type preservation
  147. A = np.array(matrix_as_list, dtype=complex)
  148. A_logm, info = logm(A, disp=False)
  149. assert_(A_logm.dtype.char in complex_dtype_chars)
  150. # check float->complex type conversion for the matrix negation
  151. A = -np.array(matrix_as_list, dtype=float)
  152. A_logm, info = logm(A, disp=False)
  153. assert_(A_logm.dtype.char in complex_dtype_chars)
  154. def test_complex_spectrum_real_logm(self):
  155. # This matrix has complex eigenvalues and real logm.
  156. # Its output dtype depends on its input dtype.
  157. M = [[1, 1, 2], [2, 1, 1], [1, 2, 1]]
  158. for dt in float, complex:
  159. X = np.array(M, dtype=dt)
  160. w = scipy.linalg.eigvals(X)
  161. assert_(1e-2 < np.absolute(w.imag).sum())
  162. Y, info = logm(X, disp=False)
  163. assert_(np.issubdtype(Y.dtype, np.inexact))
  164. assert_allclose(expm(Y), X)
  165. def test_real_mixed_sign_spectrum(self):
  166. # These matrices have real eigenvalues with mixed signs.
  167. # The output logm dtype is complex, regardless of input dtype.
  168. for M in (
  169. [[1, 0], [0, -1]],
  170. [[0, 1], [1, 0]]):
  171. for dt in float, complex:
  172. A = np.array(M, dtype=dt)
  173. A_logm, info = logm(A, disp=False)
  174. assert_(np.issubdtype(A_logm.dtype, np.complexfloating))
  175. def test_exactly_singular(self):
  176. A = np.array([[0, 0], [1j, 1j]])
  177. B = np.asarray([[1, 1], [0, 0]])
  178. for M in A, A.T, B, B.T:
  179. expected_warning = _matfuncs_inv_ssq.LogmExactlySingularWarning
  180. L, info = _assert_warns(expected_warning, logm, M, disp=False)
  181. E = expm(L)
  182. assert_allclose(E, M, atol=1e-14)
  183. def test_nearly_singular(self):
  184. M = np.array([[1e-100]])
  185. expected_warning = _matfuncs_inv_ssq.LogmNearlySingularWarning
  186. L, info = _assert_warns(expected_warning, logm, M, disp=False)
  187. E = expm(L)
  188. assert_allclose(E, M, atol=1e-14)
  189. def test_opposite_sign_complex_eigenvalues(self):
  190. # See gh-6113
  191. E = [[0, 1], [-1, 0]]
  192. L = [[0, np.pi*0.5], [-np.pi*0.5, 0]]
  193. assert_allclose(expm(L), E, atol=1e-14)
  194. assert_allclose(logm(E), L, atol=1e-14)
  195. E = [[1j, 4], [0, -1j]]
  196. L = [[1j*np.pi*0.5, 2*np.pi], [0, -1j*np.pi*0.5]]
  197. assert_allclose(expm(L), E, atol=1e-14)
  198. assert_allclose(logm(E), L, atol=1e-14)
  199. E = [[1j, 0], [0, -1j]]
  200. L = [[1j*np.pi*0.5, 0], [0, -1j*np.pi*0.5]]
  201. assert_allclose(expm(L), E, atol=1e-14)
  202. assert_allclose(logm(E), L, atol=1e-14)
  203. class TestSqrtM(object):
  204. def test_round_trip_random_float(self):
  205. np.random.seed(1234)
  206. for n in range(1, 6):
  207. M_unscaled = np.random.randn(n, n)
  208. for scale in np.logspace(-4, 4, 9):
  209. M = M_unscaled * scale
  210. M_sqrtm, info = sqrtm(M, disp=False)
  211. M_sqrtm_round_trip = M_sqrtm.dot(M_sqrtm)
  212. assert_allclose(M_sqrtm_round_trip, M)
  213. def test_round_trip_random_complex(self):
  214. np.random.seed(1234)
  215. for n in range(1, 6):
  216. M_unscaled = np.random.randn(n, n) + 1j * np.random.randn(n, n)
  217. for scale in np.logspace(-4, 4, 9):
  218. M = M_unscaled * scale
  219. M_sqrtm, info = sqrtm(M, disp=False)
  220. M_sqrtm_round_trip = M_sqrtm.dot(M_sqrtm)
  221. assert_allclose(M_sqrtm_round_trip, M)
  222. def test_bad(self):
  223. # See https://web.archive.org/web/20051220232650/http://www.maths.man.ac.uk/~nareports/narep336.ps.gz
  224. e = 2**-5
  225. se = sqrt(e)
  226. a = array([[1.0,0,0,1],
  227. [0,e,0,0],
  228. [0,0,e,0],
  229. [0,0,0,1]])
  230. sa = array([[1,0,0,0.5],
  231. [0,se,0,0],
  232. [0,0,se,0],
  233. [0,0,0,1]])
  234. n = a.shape[0]
  235. assert_array_almost_equal(dot(sa,sa),a)
  236. # Check default sqrtm.
  237. esa = sqrtm(a, disp=False, blocksize=n)[0]
  238. assert_array_almost_equal(dot(esa,esa),a)
  239. # Check sqrtm with 2x2 blocks.
  240. esa = sqrtm(a, disp=False, blocksize=2)[0]
  241. assert_array_almost_equal(dot(esa,esa),a)
  242. def test_sqrtm_type_preservation_and_conversion(self):
  243. # The sqrtm matrix function should preserve the type of a matrix
  244. # whose eigenvalues are nonnegative with zero imaginary part.
  245. # Test this preservation for variously structured matrices.
  246. complex_dtype_chars = ('F', 'D', 'G')
  247. for matrix_as_list in (
  248. [[1, 0], [0, 1]],
  249. [[1, 0], [1, 1]],
  250. [[2, 1], [1, 1]],
  251. [[2, 3], [1, 2]],
  252. [[1, 1], [1, 1]]):
  253. # check that the spectrum has the expected properties
  254. W = scipy.linalg.eigvals(matrix_as_list)
  255. assert_(not any(w.imag or w.real < 0 for w in W))
  256. # check float type preservation
  257. A = np.array(matrix_as_list, dtype=float)
  258. A_sqrtm, info = sqrtm(A, disp=False)
  259. assert_(A_sqrtm.dtype.char not in complex_dtype_chars)
  260. # check complex type preservation
  261. A = np.array(matrix_as_list, dtype=complex)
  262. A_sqrtm, info = sqrtm(A, disp=False)
  263. assert_(A_sqrtm.dtype.char in complex_dtype_chars)
  264. # check float->complex type conversion for the matrix negation
  265. A = -np.array(matrix_as_list, dtype=float)
  266. A_sqrtm, info = sqrtm(A, disp=False)
  267. assert_(A_sqrtm.dtype.char in complex_dtype_chars)
  268. def test_sqrtm_type_conversion_mixed_sign_or_complex_spectrum(self):
  269. complex_dtype_chars = ('F', 'D', 'G')
  270. for matrix_as_list in (
  271. [[1, 0], [0, -1]],
  272. [[0, 1], [1, 0]],
  273. [[0, 1, 0], [0, 0, 1], [1, 0, 0]]):
  274. # check that the spectrum has the expected properties
  275. W = scipy.linalg.eigvals(matrix_as_list)
  276. assert_(any(w.imag or w.real < 0 for w in W))
  277. # check complex->complex
  278. A = np.array(matrix_as_list, dtype=complex)
  279. A_sqrtm, info = sqrtm(A, disp=False)
  280. assert_(A_sqrtm.dtype.char in complex_dtype_chars)
  281. # check float->complex
  282. A = np.array(matrix_as_list, dtype=float)
  283. A_sqrtm, info = sqrtm(A, disp=False)
  284. assert_(A_sqrtm.dtype.char in complex_dtype_chars)
  285. def test_blocksizes(self):
  286. # Make sure I do not goof up the blocksizes when they do not divide n.
  287. np.random.seed(1234)
  288. for n in range(1, 8):
  289. A = np.random.rand(n, n) + 1j*np.random.randn(n, n)
  290. A_sqrtm_default, info = sqrtm(A, disp=False, blocksize=n)
  291. assert_allclose(A, np.linalg.matrix_power(A_sqrtm_default, 2))
  292. for blocksize in range(1, 10):
  293. A_sqrtm_new, info = sqrtm(A, disp=False, blocksize=blocksize)
  294. assert_allclose(A_sqrtm_default, A_sqrtm_new)
  295. def test_al_mohy_higham_2012_experiment_1(self):
  296. # Matrix square root of a tricky upper triangular matrix.
  297. A = _get_al_mohy_higham_2012_experiment_1()
  298. A_sqrtm, info = sqrtm(A, disp=False)
  299. A_round_trip = A_sqrtm.dot(A_sqrtm)
  300. assert_allclose(A_round_trip, A, rtol=1e-5)
  301. assert_allclose(np.tril(A_round_trip), np.tril(A))
  302. def test_strict_upper_triangular(self):
  303. # This matrix has no square root.
  304. for dt in int, float:
  305. A = np.array([
  306. [0, 3, 0, 0],
  307. [0, 0, 3, 0],
  308. [0, 0, 0, 3],
  309. [0, 0, 0, 0]], dtype=dt)
  310. A_sqrtm, info = sqrtm(A, disp=False)
  311. assert_(np.isnan(A_sqrtm).all())
  312. def test_weird_matrix(self):
  313. # The square root of matrix B exists.
  314. for dt in int, float:
  315. A = np.array([
  316. [0, 0, 1],
  317. [0, 0, 0],
  318. [0, 1, 0]], dtype=dt)
  319. B = np.array([
  320. [0, 1, 0],
  321. [0, 0, 0],
  322. [0, 0, 0]], dtype=dt)
  323. assert_array_equal(B, A.dot(A))
  324. # But scipy sqrtm is not clever enough to find it.
  325. B_sqrtm, info = sqrtm(B, disp=False)
  326. assert_(np.isnan(B_sqrtm).all())
  327. def test_disp(self):
  328. from io import StringIO
  329. np.random.seed(1234)
  330. A = np.random.rand(3, 3)
  331. B = sqrtm(A, disp=True)
  332. assert_allclose(B.dot(B), A)
  333. def test_opposite_sign_complex_eigenvalues(self):
  334. M = [[2j, 4], [0, -2j]]
  335. R = [[1+1j, 2], [0, 1-1j]]
  336. assert_allclose(np.dot(R, R), M, atol=1e-14)
  337. assert_allclose(sqrtm(M), R, atol=1e-14)
  338. def test_gh4866(self):
  339. M = np.array([[1, 0, 0, 1],
  340. [0, 0, 0, 0],
  341. [0, 0, 0, 0],
  342. [1, 0, 0, 1]])
  343. R = np.array([[sqrt(0.5), 0, 0, sqrt(0.5)],
  344. [0, 0, 0, 0],
  345. [0, 0, 0, 0],
  346. [sqrt(0.5), 0, 0, sqrt(0.5)]])
  347. assert_allclose(np.dot(R, R), M, atol=1e-14)
  348. assert_allclose(sqrtm(M), R, atol=1e-14)
  349. def test_gh5336(self):
  350. M = np.diag([2, 1, 0])
  351. R = np.diag([sqrt(2), 1, 0])
  352. assert_allclose(np.dot(R, R), M, atol=1e-14)
  353. assert_allclose(sqrtm(M), R, atol=1e-14)
  354. def test_gh7839(self):
  355. M = np.zeros((2, 2))
  356. R = np.zeros((2, 2))
  357. assert_allclose(np.dot(R, R), M, atol=1e-14)
  358. assert_allclose(sqrtm(M), R, atol=1e-14)
  359. class TestFractionalMatrixPower(object):
  360. def test_round_trip_random_complex(self):
  361. np.random.seed(1234)
  362. for p in range(1, 5):
  363. for n in range(1, 5):
  364. M_unscaled = np.random.randn(n, n) + 1j * np.random.randn(n, n)
  365. for scale in np.logspace(-4, 4, 9):
  366. M = M_unscaled * scale
  367. M_root = fractional_matrix_power(M, 1/p)
  368. M_round_trip = np.linalg.matrix_power(M_root, p)
  369. assert_allclose(M_round_trip, M)
  370. def test_round_trip_random_float(self):
  371. # This test is more annoying because it can hit the branch cut;
  372. # this happens when the matrix has an eigenvalue
  373. # with no imaginary component and with a real negative component,
  374. # and it means that the principal branch does not exist.
  375. np.random.seed(1234)
  376. for p in range(1, 5):
  377. for n in range(1, 5):
  378. M_unscaled = np.random.randn(n, n)
  379. for scale in np.logspace(-4, 4, 9):
  380. M = M_unscaled * scale
  381. M_root = fractional_matrix_power(M, 1/p)
  382. M_round_trip = np.linalg.matrix_power(M_root, p)
  383. assert_allclose(M_round_trip, M)
  384. def test_larger_abs_fractional_matrix_powers(self):
  385. np.random.seed(1234)
  386. for n in (2, 3, 5):
  387. for i in range(10):
  388. M = np.random.randn(n, n) + 1j * np.random.randn(n, n)
  389. M_one_fifth = fractional_matrix_power(M, 0.2)
  390. # Test the round trip.
  391. M_round_trip = np.linalg.matrix_power(M_one_fifth, 5)
  392. assert_allclose(M, M_round_trip)
  393. # Test a large abs fractional power.
  394. X = fractional_matrix_power(M, -5.4)
  395. Y = np.linalg.matrix_power(M_one_fifth, -27)
  396. assert_allclose(X, Y)
  397. # Test another large abs fractional power.
  398. X = fractional_matrix_power(M, 3.8)
  399. Y = np.linalg.matrix_power(M_one_fifth, 19)
  400. assert_allclose(X, Y)
  401. def test_random_matrices_and_powers(self):
  402. # Each independent iteration of this fuzz test picks random parameters.
  403. # It tries to hit some edge cases.
  404. np.random.seed(1234)
  405. nsamples = 20
  406. for i in range(nsamples):
  407. # Sample a matrix size and a random real power.
  408. n = random.randrange(1, 5)
  409. p = np.random.randn()
  410. # Sample a random real or complex matrix.
  411. matrix_scale = np.exp(random.randrange(-4, 5))
  412. A = np.random.randn(n, n)
  413. if random.choice((True, False)):
  414. A = A + 1j * np.random.randn(n, n)
  415. A = A * matrix_scale
  416. # Check a couple of analytically equivalent ways
  417. # to compute the fractional matrix power.
  418. # These can be compared because they both use the principal branch.
  419. A_power = fractional_matrix_power(A, p)
  420. A_logm, info = logm(A, disp=False)
  421. A_power_expm_logm = expm(A_logm * p)
  422. assert_allclose(A_power, A_power_expm_logm)
  423. def test_al_mohy_higham_2012_experiment_1(self):
  424. # Fractional powers of a tricky upper triangular matrix.
  425. A = _get_al_mohy_higham_2012_experiment_1()
  426. # Test remainder matrix power.
  427. A_funm_sqrt, info = funm(A, np.sqrt, disp=False)
  428. A_sqrtm, info = sqrtm(A, disp=False)
  429. A_rem_power = _matfuncs_inv_ssq._remainder_matrix_power(A, 0.5)
  430. A_power = fractional_matrix_power(A, 0.5)
  431. assert_array_equal(A_rem_power, A_power)
  432. assert_allclose(A_sqrtm, A_power)
  433. assert_allclose(A_sqrtm, A_funm_sqrt)
  434. # Test more fractional powers.
  435. for p in (1/2, 5/3):
  436. A_power = fractional_matrix_power(A, p)
  437. A_round_trip = fractional_matrix_power(A_power, 1/p)
  438. assert_allclose(A_round_trip, A, rtol=1e-2)
  439. assert_allclose(np.tril(A_round_trip, 1), np.tril(A, 1))
  440. def test_briggs_helper_function(self):
  441. np.random.seed(1234)
  442. for a in np.random.randn(10) + 1j * np.random.randn(10):
  443. for k in range(5):
  444. x_observed = _matfuncs_inv_ssq._briggs_helper_function(a, k)
  445. x_expected = a ** np.exp2(-k) - 1
  446. assert_allclose(x_observed, x_expected)
  447. def test_type_preservation_and_conversion(self):
  448. # The fractional_matrix_power matrix function should preserve
  449. # the type of a matrix whose eigenvalues
  450. # are positive with zero imaginary part.
  451. # Test this preservation for variously structured matrices.
  452. complex_dtype_chars = ('F', 'D', 'G')
  453. for matrix_as_list in (
  454. [[1, 0], [0, 1]],
  455. [[1, 0], [1, 1]],
  456. [[2, 1], [1, 1]],
  457. [[2, 3], [1, 2]]):
  458. # check that the spectrum has the expected properties
  459. W = scipy.linalg.eigvals(matrix_as_list)
  460. assert_(not any(w.imag or w.real < 0 for w in W))
  461. # Check various positive and negative powers
  462. # with absolute values bigger and smaller than 1.
  463. for p in (-2.4, -0.9, 0.2, 3.3):
  464. # check float type preservation
  465. A = np.array(matrix_as_list, dtype=float)
  466. A_power = fractional_matrix_power(A, p)
  467. assert_(A_power.dtype.char not in complex_dtype_chars)
  468. # check complex type preservation
  469. A = np.array(matrix_as_list, dtype=complex)
  470. A_power = fractional_matrix_power(A, p)
  471. assert_(A_power.dtype.char in complex_dtype_chars)
  472. # check float->complex for the matrix negation
  473. A = -np.array(matrix_as_list, dtype=float)
  474. A_power = fractional_matrix_power(A, p)
  475. assert_(A_power.dtype.char in complex_dtype_chars)
  476. def test_type_conversion_mixed_sign_or_complex_spectrum(self):
  477. complex_dtype_chars = ('F', 'D', 'G')
  478. for matrix_as_list in (
  479. [[1, 0], [0, -1]],
  480. [[0, 1], [1, 0]],
  481. [[0, 1, 0], [0, 0, 1], [1, 0, 0]]):
  482. # check that the spectrum has the expected properties
  483. W = scipy.linalg.eigvals(matrix_as_list)
  484. assert_(any(w.imag or w.real < 0 for w in W))
  485. # Check various positive and negative powers
  486. # with absolute values bigger and smaller than 1.
  487. for p in (-2.4, -0.9, 0.2, 3.3):
  488. # check complex->complex
  489. A = np.array(matrix_as_list, dtype=complex)
  490. A_power = fractional_matrix_power(A, p)
  491. assert_(A_power.dtype.char in complex_dtype_chars)
  492. # check float->complex
  493. A = np.array(matrix_as_list, dtype=float)
  494. A_power = fractional_matrix_power(A, p)
  495. assert_(A_power.dtype.char in complex_dtype_chars)
  496. @pytest.mark.xfail(reason='Too unstable across LAPACKs.')
  497. def test_singular(self):
  498. # Negative fractional powers do not work with singular matrices.
  499. for matrix_as_list in (
  500. [[0, 0], [0, 0]],
  501. [[1, 1], [1, 1]],
  502. [[1, 2], [3, 6]],
  503. [[0, 0, 0], [0, 1, 1], [0, -1, 1]]):
  504. # Check fractional powers both for float and for complex types.
  505. for newtype in (float, complex):
  506. A = np.array(matrix_as_list, dtype=newtype)
  507. for p in (-0.7, -0.9, -2.4, -1.3):
  508. A_power = fractional_matrix_power(A, p)
  509. assert_(np.isnan(A_power).all())
  510. for p in (0.2, 1.43):
  511. A_power = fractional_matrix_power(A, p)
  512. A_round_trip = fractional_matrix_power(A_power, 1/p)
  513. assert_allclose(A_round_trip, A)
  514. def test_opposite_sign_complex_eigenvalues(self):
  515. M = [[2j, 4], [0, -2j]]
  516. R = [[1+1j, 2], [0, 1-1j]]
  517. assert_allclose(np.dot(R, R), M, atol=1e-14)
  518. assert_allclose(fractional_matrix_power(M, 0.5), R, atol=1e-14)
  519. class TestExpM(object):
  520. def test_zero(self):
  521. a = array([[0.,0],[0,0]])
  522. assert_array_almost_equal(expm(a),[[1,0],[0,1]])
  523. def test_single_elt(self):
  524. # See gh-5853
  525. from scipy.sparse import csc_matrix
  526. vOne = -2.02683397006j
  527. vTwo = -2.12817566856j
  528. mOne = csc_matrix([[vOne]], dtype='complex')
  529. mTwo = csc_matrix([[vTwo]], dtype='complex')
  530. outOne = expm(mOne)
  531. outTwo = expm(mTwo)
  532. assert_equal(type(outOne), type(mOne))
  533. assert_equal(type(outTwo), type(mTwo))
  534. assert_allclose(outOne[0, 0], complex(-0.44039415155949196,
  535. -0.8978045395698304))
  536. assert_allclose(outTwo[0, 0], complex(-0.52896401032626006,
  537. -0.84864425749518878))
  538. class TestExpmFrechet(object):
  539. def test_expm_frechet(self):
  540. # a test of the basic functionality
  541. M = np.array([
  542. [1, 2, 3, 4],
  543. [5, 6, 7, 8],
  544. [0, 0, 1, 2],
  545. [0, 0, 5, 6],
  546. ], dtype=float)
  547. A = np.array([
  548. [1, 2],
  549. [5, 6],
  550. ], dtype=float)
  551. E = np.array([
  552. [3, 4],
  553. [7, 8],
  554. ], dtype=float)
  555. expected_expm = scipy.linalg.expm(A)
  556. expected_frechet = scipy.linalg.expm(M)[:2, 2:]
  557. for kwargs in ({}, {'method':'SPS'}, {'method':'blockEnlarge'}):
  558. observed_expm, observed_frechet = expm_frechet(A, E, **kwargs)
  559. assert_allclose(expected_expm, observed_expm)
  560. assert_allclose(expected_frechet, observed_frechet)
  561. def test_small_norm_expm_frechet(self):
  562. # methodically test matrices with a range of norms, for better coverage
  563. M_original = np.array([
  564. [1, 2, 3, 4],
  565. [5, 6, 7, 8],
  566. [0, 0, 1, 2],
  567. [0, 0, 5, 6],
  568. ], dtype=float)
  569. A_original = np.array([
  570. [1, 2],
  571. [5, 6],
  572. ], dtype=float)
  573. E_original = np.array([
  574. [3, 4],
  575. [7, 8],
  576. ], dtype=float)
  577. A_original_norm_1 = scipy.linalg.norm(A_original, 1)
  578. selected_m_list = [1, 3, 5, 7, 9, 11, 13, 15]
  579. m_neighbor_pairs = zip(selected_m_list[:-1], selected_m_list[1:])
  580. for ma, mb in m_neighbor_pairs:
  581. ell_a = scipy.linalg._expm_frechet.ell_table_61[ma]
  582. ell_b = scipy.linalg._expm_frechet.ell_table_61[mb]
  583. target_norm_1 = 0.5 * (ell_a + ell_b)
  584. scale = target_norm_1 / A_original_norm_1
  585. M = scale * M_original
  586. A = scale * A_original
  587. E = scale * E_original
  588. expected_expm = scipy.linalg.expm(A)
  589. expected_frechet = scipy.linalg.expm(M)[:2, 2:]
  590. observed_expm, observed_frechet = expm_frechet(A, E)
  591. assert_allclose(expected_expm, observed_expm)
  592. assert_allclose(expected_frechet, observed_frechet)
  593. def test_fuzz(self):
  594. # try a bunch of crazy inputs
  595. rfuncs = (
  596. np.random.uniform,
  597. np.random.normal,
  598. np.random.standard_cauchy,
  599. np.random.exponential)
  600. ntests = 100
  601. for i in range(ntests):
  602. rfunc = random.choice(rfuncs)
  603. target_norm_1 = random.expovariate(1.0)
  604. n = random.randrange(2, 16)
  605. A_original = rfunc(size=(n,n))
  606. E_original = rfunc(size=(n,n))
  607. A_original_norm_1 = scipy.linalg.norm(A_original, 1)
  608. scale = target_norm_1 / A_original_norm_1
  609. A = scale * A_original
  610. E = scale * E_original
  611. M = np.vstack([
  612. np.hstack([A, E]),
  613. np.hstack([np.zeros_like(A), A])])
  614. expected_expm = scipy.linalg.expm(A)
  615. expected_frechet = scipy.linalg.expm(M)[:n, n:]
  616. observed_expm, observed_frechet = expm_frechet(A, E)
  617. assert_allclose(expected_expm, observed_expm)
  618. assert_allclose(expected_frechet, observed_frechet)
  619. def test_problematic_matrix(self):
  620. # this test case uncovered a bug which has since been fixed
  621. A = np.array([
  622. [1.50591997, 1.93537998],
  623. [0.41203263, 0.23443516],
  624. ], dtype=float)
  625. E = np.array([
  626. [1.87864034, 2.07055038],
  627. [1.34102727, 0.67341123],
  628. ], dtype=float)
  629. A_norm_1 = scipy.linalg.norm(A, 1)
  630. sps_expm, sps_frechet = expm_frechet(
  631. A, E, method='SPS')
  632. blockEnlarge_expm, blockEnlarge_frechet = expm_frechet(
  633. A, E, method='blockEnlarge')
  634. assert_allclose(sps_expm, blockEnlarge_expm)
  635. assert_allclose(sps_frechet, blockEnlarge_frechet)
  636. @pytest.mark.slow
  637. @pytest.mark.skip(reason='this test is deliberately slow')
  638. def test_medium_matrix(self):
  639. # profile this to see the speed difference
  640. n = 1000
  641. A = np.random.exponential(size=(n, n))
  642. E = np.random.exponential(size=(n, n))
  643. sps_expm, sps_frechet = expm_frechet(
  644. A, E, method='SPS')
  645. blockEnlarge_expm, blockEnlarge_frechet = expm_frechet(
  646. A, E, method='blockEnlarge')
  647. assert_allclose(sps_expm, blockEnlarge_expm)
  648. assert_allclose(sps_frechet, blockEnlarge_frechet)
  649. def _help_expm_cond_search(A, A_norm, X, X_norm, eps, p):
  650. p = np.reshape(p, A.shape)
  651. p_norm = norm(p)
  652. perturbation = eps * p * (A_norm / p_norm)
  653. X_prime = expm(A + perturbation)
  654. scaled_relative_error = norm(X_prime - X) / (X_norm * eps)
  655. return -scaled_relative_error
  656. def _normalized_like(A, B):
  657. return A * (scipy.linalg.norm(B) / scipy.linalg.norm(A))
  658. def _relative_error(f, A, perturbation):
  659. X = f(A)
  660. X_prime = f(A + perturbation)
  661. return norm(X_prime - X) / norm(X)
  662. class TestExpmConditionNumber(object):
  663. def test_expm_cond_smoke(self):
  664. np.random.seed(1234)
  665. for n in range(1, 4):
  666. A = np.random.randn(n, n)
  667. kappa = expm_cond(A)
  668. assert_array_less(0, kappa)
  669. def test_expm_bad_condition_number(self):
  670. A = np.array([
  671. [-1.128679820, 9.614183771e4, -4.524855739e9, 2.924969411e14],
  672. [0, -1.201010529, 9.634696872e4, -4.681048289e9],
  673. [0, 0, -1.132893222, 9.532491830e4],
  674. [0, 0, 0, -1.179475332],
  675. ])
  676. kappa = expm_cond(A)
  677. assert_array_less(1e36, kappa)
  678. def test_univariate(self):
  679. np.random.seed(12345)
  680. for x in np.linspace(-5, 5, num=11):
  681. A = np.array([[x]])
  682. assert_allclose(expm_cond(A), abs(x))
  683. for x in np.logspace(-2, 2, num=11):
  684. A = np.array([[x]])
  685. assert_allclose(expm_cond(A), abs(x))
  686. for i in range(10):
  687. A = np.random.randn(1, 1)
  688. assert_allclose(expm_cond(A), np.absolute(A)[0, 0])
  689. @pytest.mark.slow
  690. def test_expm_cond_fuzz(self):
  691. np.random.seed(12345)
  692. eps = 1e-5
  693. nsamples = 10
  694. for i in range(nsamples):
  695. n = np.random.randint(2, 5)
  696. A = np.random.randn(n, n)
  697. A_norm = scipy.linalg.norm(A)
  698. X = expm(A)
  699. X_norm = scipy.linalg.norm(X)
  700. kappa = expm_cond(A)
  701. # Look for the small perturbation that gives the greatest
  702. # relative error.
  703. f = functools.partial(_help_expm_cond_search,
  704. A, A_norm, X, X_norm, eps)
  705. guess = np.ones(n*n)
  706. out = minimize(f, guess, method='L-BFGS-B')
  707. xopt = out.x
  708. yopt = f(xopt)
  709. p_best = eps * _normalized_like(np.reshape(xopt, A.shape), A)
  710. p_best_relerr = _relative_error(expm, A, p_best)
  711. assert_allclose(p_best_relerr, -yopt * eps)
  712. # Check that the identified perturbation indeed gives greater
  713. # relative error than random perturbations with similar norms.
  714. for j in range(5):
  715. p_rand = eps * _normalized_like(np.random.randn(*A.shape), A)
  716. assert_allclose(norm(p_best), norm(p_rand))
  717. p_rand_relerr = _relative_error(expm, A, p_rand)
  718. assert_array_less(p_rand_relerr, p_best_relerr)
  719. # The greatest relative error should not be much greater than
  720. # eps times the condition number kappa.
  721. # In the limit as eps approaches zero it should never be greater.
  722. assert_array_less(p_best_relerr, (1 + 2*eps) * eps * kappa)