test_regression.py 6.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190
  1. """ Test functions for linalg module
  2. """
  3. import pytest
  4. import numpy as np
  5. from numpy import arange, array, dot, float64, linalg, transpose
  6. from numpy.testing import (
  7. assert_,
  8. assert_almost_equal,
  9. assert_array_almost_equal,
  10. assert_array_equal,
  11. assert_array_less,
  12. assert_equal,
  13. assert_raises,
  14. )
  15. class TestRegression:
  16. def test_eig_build(self):
  17. # Ticket #652
  18. rva = array([1.03221168e+02 + 0.j,
  19. -1.91843603e+01 + 0.j,
  20. -6.04004526e-01 + 15.84422474j,
  21. -6.04004526e-01 - 15.84422474j,
  22. -1.13692929e+01 + 0.j,
  23. -6.57612485e-01 + 10.41755503j,
  24. -6.57612485e-01 - 10.41755503j,
  25. 1.82126812e+01 + 0.j,
  26. 1.06011014e+01 + 0.j,
  27. 7.80732773e+00 + 0.j,
  28. -7.65390898e-01 + 0.j,
  29. 1.51971555e-15 + 0.j,
  30. -1.51308713e-15 + 0.j])
  31. a = arange(13 * 13, dtype=float64)
  32. a = a.reshape((13, 13))
  33. a = a % 17
  34. va, ve = linalg.eig(a)
  35. va.sort()
  36. rva.sort()
  37. assert_array_almost_equal(va, rva)
  38. def test_eigh_build(self):
  39. # Ticket 662.
  40. rvals = [68.60568999, 89.57756725, 106.67185574]
  41. cov = array([[77.70273908, 3.51489954, 15.64602427],
  42. [ 3.51489954, 88.97013878, -1.07431931],
  43. [15.64602427, -1.07431931, 98.18223512]])
  44. vals, vecs = linalg.eigh(cov)
  45. assert_array_almost_equal(vals, rvals)
  46. def test_svd_build(self):
  47. # Ticket 627.
  48. a = array([[0., 1.], [1., 1.], [2., 1.], [3., 1.]])
  49. m, n = a.shape
  50. u, s, vh = linalg.svd(a)
  51. b = dot(transpose(u[:, n:]), a)
  52. assert_array_almost_equal(b, np.zeros((2, 2)))
  53. def test_norm_vector_badarg(self):
  54. # Regression for #786: Frobenius norm for vectors raises
  55. # ValueError.
  56. assert_raises(ValueError, linalg.norm, array([1., 2., 3.]), 'fro')
  57. def test_lapack_endian(self):
  58. # For bug #1482
  59. a = array([[ 5.7998084, -2.1825367],
  60. [-2.1825367, 9.85910595]], dtype='>f8')
  61. b = array(a, dtype='<f8')
  62. ap = linalg.cholesky(a)
  63. bp = linalg.cholesky(b)
  64. assert_array_equal(ap, bp)
  65. def test_large_svd_32bit(self):
  66. # See gh-4442, 64bit would require very large/slow matrices.
  67. x = np.eye(1000, 66)
  68. np.linalg.svd(x)
  69. def test_svd_no_uv(self):
  70. # gh-4733
  71. for shape in (3, 4), (4, 4), (4, 3):
  72. for t in float, complex:
  73. a = np.ones(shape, dtype=t)
  74. w = linalg.svd(a, compute_uv=False)
  75. c = np.count_nonzero(np.absolute(w) > 0.5)
  76. assert_equal(c, 1)
  77. assert_equal(np.linalg.matrix_rank(a), 1)
  78. assert_array_less(1, np.linalg.norm(a, ord=2))
  79. w_svdvals = linalg.svdvals(a)
  80. assert_array_almost_equal(w, w_svdvals)
  81. def test_norm_object_array(self):
  82. # gh-7575
  83. testvector = np.array([np.array([0, 1]), 0, 0], dtype=object)
  84. norm = linalg.norm(testvector)
  85. assert_array_equal(norm, [0, 1])
  86. assert_(norm.dtype == np.dtype('float64'))
  87. norm = linalg.norm(testvector, ord=1)
  88. assert_array_equal(norm, [0, 1])
  89. assert_(norm.dtype != np.dtype('float64'))
  90. norm = linalg.norm(testvector, ord=2)
  91. assert_array_equal(norm, [0, 1])
  92. assert_(norm.dtype == np.dtype('float64'))
  93. assert_raises(ValueError, linalg.norm, testvector, ord='fro')
  94. assert_raises(ValueError, linalg.norm, testvector, ord='nuc')
  95. assert_raises(ValueError, linalg.norm, testvector, ord=np.inf)
  96. assert_raises(ValueError, linalg.norm, testvector, ord=-np.inf)
  97. assert_raises(ValueError, linalg.norm, testvector, ord=0)
  98. assert_raises(ValueError, linalg.norm, testvector, ord=-1)
  99. assert_raises(ValueError, linalg.norm, testvector, ord=-2)
  100. testmatrix = np.array([[np.array([0, 1]), 0, 0],
  101. [0, 0, 0]], dtype=object)
  102. norm = linalg.norm(testmatrix)
  103. assert_array_equal(norm, [0, 1])
  104. assert_(norm.dtype == np.dtype('float64'))
  105. norm = linalg.norm(testmatrix, ord='fro')
  106. assert_array_equal(norm, [0, 1])
  107. assert_(norm.dtype == np.dtype('float64'))
  108. assert_raises(TypeError, linalg.norm, testmatrix, ord='nuc')
  109. assert_raises(ValueError, linalg.norm, testmatrix, ord=np.inf)
  110. assert_raises(ValueError, linalg.norm, testmatrix, ord=-np.inf)
  111. assert_raises(ValueError, linalg.norm, testmatrix, ord=0)
  112. assert_raises(ValueError, linalg.norm, testmatrix, ord=1)
  113. assert_raises(ValueError, linalg.norm, testmatrix, ord=-1)
  114. assert_raises(TypeError, linalg.norm, testmatrix, ord=2)
  115. assert_raises(TypeError, linalg.norm, testmatrix, ord=-2)
  116. assert_raises(ValueError, linalg.norm, testmatrix, ord=3)
  117. def test_lstsq_complex_larger_rhs(self):
  118. # gh-9891
  119. size = 20
  120. n_rhs = 70
  121. G = np.random.randn(size, size) + 1j * np.random.randn(size, size)
  122. u = np.random.randn(size, n_rhs) + 1j * np.random.randn(size, n_rhs)
  123. b = G.dot(u)
  124. # This should work without segmentation fault.
  125. u_lstsq, res, rank, sv = linalg.lstsq(G, b, rcond=None)
  126. # check results just in case
  127. assert_array_almost_equal(u_lstsq, u)
  128. @pytest.mark.parametrize("upper", [True, False])
  129. def test_cholesky_empty_array(self, upper):
  130. # gh-25840 - upper=True hung before.
  131. res = np.linalg.cholesky(np.zeros((0, 0)), upper=upper)
  132. assert res.size == 0
  133. @pytest.mark.parametrize("rtol", [0.0, [0.0] * 4, np.zeros((4,))])
  134. def test_matrix_rank_rtol_argument(self, rtol):
  135. # gh-25877
  136. x = np.zeros((4, 3, 2))
  137. res = np.linalg.matrix_rank(x, rtol=rtol)
  138. assert res.shape == (4,)
  139. @pytest.mark.thread_unsafe(reason="test is already testing threads with openblas")
  140. def test_openblas_threading(self):
  141. # gh-27036
  142. # Test whether matrix multiplication involving a large matrix always
  143. # gives the same (correct) answer
  144. x = np.arange(500000, dtype=np.float64)
  145. src = np.vstack((x, -10 * x)).T
  146. matrix = np.array([[0, 1], [1, 0]])
  147. expected = np.vstack((-10 * x, x)).T # src @ matrix
  148. for i in range(200):
  149. result = src @ matrix
  150. mismatches = (~np.isclose(result, expected)).sum()
  151. if mismatches != 0:
  152. assert False, ("unexpected result from matmul, "
  153. "probably due to OpenBLAS threading issues")
  154. def test_norm_linux_arm(self):
  155. # gh-30816
  156. a = np.arange(20000) / 50000
  157. b = a + 1j * np.roll(np.flip(a), 12345)
  158. norm = np.linalg.norm(b)
  159. assert_almost_equal(norm, 46.18628948075393)