test_matfuncs.py 38 KB

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