test_trustregion_exact.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352
  1. """
  2. Unit tests for trust-region iterative subproblem.
  3. To run it in its simplest form::
  4. nosetests test_optimize.py
  5. """
  6. import numpy as np
  7. from scipy.optimize._trustregion_exact import (
  8. estimate_smallest_singular_value,
  9. singular_leading_submatrix,
  10. IterativeSubproblem)
  11. from scipy.linalg import (svd, get_lapack_funcs, det, qr, norm)
  12. from numpy.testing import (assert_array_equal,
  13. assert_equal, assert_array_almost_equal)
  14. def random_entry(n, min_eig, max_eig, case):
  15. # Generate random matrix
  16. rand = np.random.uniform(-1, 1, (n, n))
  17. # QR decomposition
  18. Q, _, _ = qr(rand, pivoting='True')
  19. # Generate random eigenvalues
  20. eigvalues = np.random.uniform(min_eig, max_eig, n)
  21. eigvalues = np.sort(eigvalues)[::-1]
  22. # Generate matrix
  23. Qaux = np.multiply(eigvalues, Q)
  24. A = np.dot(Qaux, Q.T)
  25. # Generate gradient vector accordingly
  26. # to the case is being tested.
  27. if case == 'hard':
  28. g = np.zeros(n)
  29. g[:-1] = np.random.uniform(-1, 1, n-1)
  30. g = np.dot(Q, g)
  31. elif case == 'jac_equal_zero':
  32. g = np.zeros(n)
  33. else:
  34. g = np.random.uniform(-1, 1, n)
  35. return A, g
  36. class TestEstimateSmallestSingularValue:
  37. def test_for_ill_condiotioned_matrix(self):
  38. # Ill-conditioned triangular matrix
  39. C = np.array([[1, 2, 3, 4],
  40. [0, 0.05, 60, 7],
  41. [0, 0, 0.8, 9],
  42. [0, 0, 0, 10]])
  43. # Get svd decomposition
  44. U, s, Vt = svd(C)
  45. # Get smallest singular value and correspondent right singular vector.
  46. smin_svd = s[-1]
  47. zmin_svd = Vt[-1, :]
  48. # Estimate smallest singular value
  49. smin, zmin = estimate_smallest_singular_value(C)
  50. # Check the estimation
  51. assert_array_almost_equal(smin, smin_svd, decimal=8)
  52. assert_array_almost_equal(abs(zmin), abs(zmin_svd), decimal=8)
  53. class TestSingularLeadingSubmatrix:
  54. def test_for_already_singular_leading_submatrix(self):
  55. # Define test matrix A.
  56. # Note that the leading 2x2 submatrix is singular.
  57. A = np.array([[1, 2, 3],
  58. [2, 4, 5],
  59. [3, 5, 6]])
  60. # Get Cholesky from lapack functions
  61. cholesky, = get_lapack_funcs(('potrf',), (A,))
  62. # Compute Cholesky Decomposition
  63. c, k = cholesky(A, lower=False, overwrite_a=False, clean=True)
  64. delta, v = singular_leading_submatrix(A, c, k)
  65. A[k-1, k-1] += delta
  66. # Check if the leading submatrix is singular.
  67. assert_array_almost_equal(det(A[:k, :k]), 0)
  68. # Check if `v` fullfil the specified properties
  69. quadratic_term = np.dot(v, np.dot(A, v))
  70. assert_array_almost_equal(quadratic_term, 0)
  71. def test_for_simetric_indefinite_matrix(self):
  72. # Define test matrix A.
  73. # Note that the leading 5x5 submatrix is indefinite.
  74. A = np.asarray([[1, 2, 3, 7, 8],
  75. [2, 5, 5, 9, 0],
  76. [3, 5, 11, 1, 2],
  77. [7, 9, 1, 7, 5],
  78. [8, 0, 2, 5, 8]])
  79. # Get Cholesky from lapack functions
  80. cholesky, = get_lapack_funcs(('potrf',), (A,))
  81. # Compute Cholesky Decomposition
  82. c, k = cholesky(A, lower=False, overwrite_a=False, clean=True)
  83. delta, v = singular_leading_submatrix(A, c, k)
  84. A[k-1, k-1] += delta
  85. # Check if the leading submatrix is singular.
  86. assert_array_almost_equal(det(A[:k, :k]), 0)
  87. # Check if `v` fullfil the specified properties
  88. quadratic_term = np.dot(v, np.dot(A, v))
  89. assert_array_almost_equal(quadratic_term, 0)
  90. def test_for_first_element_equal_to_zero(self):
  91. # Define test matrix A.
  92. # Note that the leading 2x2 submatrix is singular.
  93. A = np.array([[0, 3, 11],
  94. [3, 12, 5],
  95. [11, 5, 6]])
  96. # Get Cholesky from lapack functions
  97. cholesky, = get_lapack_funcs(('potrf',), (A,))
  98. # Compute Cholesky Decomposition
  99. c, k = cholesky(A, lower=False, overwrite_a=False, clean=True)
  100. delta, v = singular_leading_submatrix(A, c, k)
  101. A[k-1, k-1] += delta
  102. # Check if the leading submatrix is singular
  103. assert_array_almost_equal(det(A[:k, :k]), 0)
  104. # Check if `v` fullfil the specified properties
  105. quadratic_term = np.dot(v, np.dot(A, v))
  106. assert_array_almost_equal(quadratic_term, 0)
  107. class TestIterativeSubproblem:
  108. def test_for_the_easy_case(self):
  109. # `H` is chosen such that `g` is not orthogonal to the
  110. # eigenvector associated with the smallest eigenvalue `s`.
  111. H = [[10, 2, 3, 4],
  112. [2, 1, 7, 1],
  113. [3, 7, 1, 7],
  114. [4, 1, 7, 2]]
  115. g = [1, 1, 1, 1]
  116. # Trust Radius
  117. trust_radius = 1
  118. # Solve Subproblem
  119. subprob = IterativeSubproblem(x=0,
  120. fun=lambda x: 0,
  121. jac=lambda x: np.array(g),
  122. hess=lambda x: np.array(H),
  123. k_easy=1e-10,
  124. k_hard=1e-10)
  125. p, hits_boundary = subprob.solve(trust_radius)
  126. assert_array_almost_equal(p, [0.00393332, -0.55260862,
  127. 0.67065477, -0.49480341])
  128. assert_array_almost_equal(hits_boundary, True)
  129. def test_for_the_hard_case(self):
  130. # `H` is chosen such that `g` is orthogonal to the
  131. # eigenvector associated with the smallest eigenvalue `s`.
  132. H = [[10, 2, 3, 4],
  133. [2, 1, 7, 1],
  134. [3, 7, 1, 7],
  135. [4, 1, 7, 2]]
  136. g = [6.4852641521327437, 1, 1, 1]
  137. s = -8.2151519874416614
  138. # Trust Radius
  139. trust_radius = 1
  140. # Solve Subproblem
  141. subprob = IterativeSubproblem(x=0,
  142. fun=lambda x: 0,
  143. jac=lambda x: np.array(g),
  144. hess=lambda x: np.array(H),
  145. k_easy=1e-10,
  146. k_hard=1e-10)
  147. p, hits_boundary = subprob.solve(trust_radius)
  148. assert_array_almost_equal(-s, subprob.lambda_current)
  149. def test_for_interior_convergence(self):
  150. H = [[1.812159, 0.82687265, 0.21838879, -0.52487006, 0.25436988],
  151. [0.82687265, 2.66380283, 0.31508988, -0.40144163, 0.08811588],
  152. [0.21838879, 0.31508988, 2.38020726, -0.3166346, 0.27363867],
  153. [-0.52487006, -0.40144163, -0.3166346, 1.61927182, -0.42140166],
  154. [0.25436988, 0.08811588, 0.27363867, -0.42140166, 1.33243101]]
  155. g = [0.75798952, 0.01421945, 0.33847612, 0.83725004, -0.47909534]
  156. # Solve Subproblem
  157. subprob = IterativeSubproblem(x=0,
  158. fun=lambda x: 0,
  159. jac=lambda x: np.array(g),
  160. hess=lambda x: np.array(H))
  161. p, hits_boundary = subprob.solve(1.1)
  162. assert_array_almost_equal(p, [-0.68585435, 0.1222621, -0.22090999,
  163. -0.67005053, 0.31586769])
  164. assert_array_almost_equal(hits_boundary, False)
  165. assert_array_almost_equal(subprob.lambda_current, 0)
  166. assert_array_almost_equal(subprob.niter, 1)
  167. def test_for_jac_equal_zero(self):
  168. H = [[0.88547534, 2.90692271, 0.98440885, -0.78911503, -0.28035809],
  169. [2.90692271, -0.04618819, 0.32867263, -0.83737945, 0.17116396],
  170. [0.98440885, 0.32867263, -0.87355957, -0.06521957, -1.43030957],
  171. [-0.78911503, -0.83737945, -0.06521957, -1.645709, -0.33887298],
  172. [-0.28035809, 0.17116396, -1.43030957, -0.33887298, -1.68586978]]
  173. g = [0, 0, 0, 0, 0]
  174. # Solve Subproblem
  175. subprob = IterativeSubproblem(x=0,
  176. fun=lambda x: 0,
  177. jac=lambda x: np.array(g),
  178. hess=lambda x: np.array(H),
  179. k_easy=1e-10,
  180. k_hard=1e-10)
  181. p, hits_boundary = subprob.solve(1.1)
  182. assert_array_almost_equal(p, [0.06910534, -0.01432721,
  183. -0.65311947, -0.23815972,
  184. -0.84954934])
  185. assert_array_almost_equal(hits_boundary, True)
  186. def test_for_jac_very_close_to_zero(self):
  187. H = [[0.88547534, 2.90692271, 0.98440885, -0.78911503, -0.28035809],
  188. [2.90692271, -0.04618819, 0.32867263, -0.83737945, 0.17116396],
  189. [0.98440885, 0.32867263, -0.87355957, -0.06521957, -1.43030957],
  190. [-0.78911503, -0.83737945, -0.06521957, -1.645709, -0.33887298],
  191. [-0.28035809, 0.17116396, -1.43030957, -0.33887298, -1.68586978]]
  192. g = [0, 0, 0, 0, 1e-15]
  193. # Solve Subproblem
  194. subprob = IterativeSubproblem(x=0,
  195. fun=lambda x: 0,
  196. jac=lambda x: np.array(g),
  197. hess=lambda x: np.array(H),
  198. k_easy=1e-10,
  199. k_hard=1e-10)
  200. p, hits_boundary = subprob.solve(1.1)
  201. assert_array_almost_equal(p, [0.06910534, -0.01432721,
  202. -0.65311947, -0.23815972,
  203. -0.84954934])
  204. assert_array_almost_equal(hits_boundary, True)
  205. def test_for_random_entries(self):
  206. # Seed
  207. np.random.seed(1)
  208. # Dimension
  209. n = 5
  210. for case in ('easy', 'hard', 'jac_equal_zero'):
  211. eig_limits = [(-20, -15),
  212. (-10, -5),
  213. (-10, 0),
  214. (-5, 5),
  215. (-10, 10),
  216. (0, 10),
  217. (5, 10),
  218. (15, 20)]
  219. for min_eig, max_eig in eig_limits:
  220. # Generate random symmetric matrix H with
  221. # eigenvalues between min_eig and max_eig.
  222. H, g = random_entry(n, min_eig, max_eig, case)
  223. # Trust radius
  224. trust_radius_list = [0.1, 0.3, 0.6, 0.8, 1, 1.2, 3.3, 5.5, 10]
  225. for trust_radius in trust_radius_list:
  226. # Solve subproblem with very high accuracy
  227. subprob_ac = IterativeSubproblem(0,
  228. lambda x: 0,
  229. lambda x: g,
  230. lambda x: H,
  231. k_easy=1e-10,
  232. k_hard=1e-10)
  233. p_ac, hits_boundary_ac = subprob_ac.solve(trust_radius)
  234. # Compute objective function value
  235. J_ac = 1/2*np.dot(p_ac, np.dot(H, p_ac))+np.dot(g, p_ac)
  236. stop_criteria = [(0.1, 2),
  237. (0.5, 1.1),
  238. (0.9, 1.01)]
  239. for k_opt, k_trf in stop_criteria:
  240. # k_easy and k_hard computed in function
  241. # of k_opt and k_trf accordingly to
  242. # Conn, A. R., Gould, N. I., & Toint, P. L. (2000).
  243. # "Trust region methods". Siam. p. 197.
  244. k_easy = min(k_trf-1,
  245. 1-np.sqrt(k_opt))
  246. k_hard = 1-k_opt
  247. # Solve subproblem
  248. subprob = IterativeSubproblem(0,
  249. lambda x: 0,
  250. lambda x: g,
  251. lambda x: H,
  252. k_easy=k_easy,
  253. k_hard=k_hard)
  254. p, hits_boundary = subprob.solve(trust_radius)
  255. # Compute objective function value
  256. J = 1/2*np.dot(p, np.dot(H, p))+np.dot(g, p)
  257. # Check if it respect k_trf
  258. if hits_boundary:
  259. assert_array_equal(np.abs(norm(p)-trust_radius) <=
  260. (k_trf-1)*trust_radius, True)
  261. else:
  262. assert_equal(norm(p) <= trust_radius, True)
  263. # Check if it respect k_opt
  264. assert_equal(J <= k_opt*J_ac, True)