test_trustregion_exact.py 13 KB

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