test_nonlin.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448
  1. """ Unit tests for nonlinear solvers
  2. Author: Ondrej Certik
  3. May 2007
  4. """
  5. from __future__ import division, print_function, absolute_import
  6. from numpy.testing import assert_
  7. import pytest
  8. from scipy._lib.six import xrange
  9. from scipy.optimize import nonlin, root
  10. from numpy import matrix, diag, dot
  11. from numpy.linalg import inv
  12. import numpy as np
  13. from .test_minpack import pressure_network
  14. SOLVERS = {'anderson': nonlin.anderson, 'diagbroyden': nonlin.diagbroyden,
  15. 'linearmixing': nonlin.linearmixing, 'excitingmixing': nonlin.excitingmixing,
  16. 'broyden1': nonlin.broyden1, 'broyden2': nonlin.broyden2,
  17. 'krylov': nonlin.newton_krylov}
  18. MUST_WORK = {'anderson': nonlin.anderson, 'broyden1': nonlin.broyden1,
  19. 'broyden2': nonlin.broyden2, 'krylov': nonlin.newton_krylov}
  20. #-------------------------------------------------------------------------------
  21. # Test problems
  22. #-------------------------------------------------------------------------------
  23. def F(x):
  24. x = np.asmatrix(x).T
  25. d = matrix(diag([3,2,1.5,1,0.5]))
  26. c = 0.01
  27. f = -d*x - c*float(x.T*x)*x
  28. return f
  29. F.xin = [1,1,1,1,1]
  30. F.KNOWN_BAD = {}
  31. def F2(x):
  32. return x
  33. F2.xin = [1,2,3,4,5,6]
  34. F2.KNOWN_BAD = {'linearmixing': nonlin.linearmixing,
  35. 'excitingmixing': nonlin.excitingmixing}
  36. def F2_lucky(x):
  37. return x
  38. F2_lucky.xin = [0,0,0,0,0,0]
  39. F2_lucky.KNOWN_BAD = {}
  40. def F3(x):
  41. A = np.mat('-2 1 0; 1 -2 1; 0 1 -2')
  42. b = np.mat('1 2 3')
  43. return np.dot(A, x) - b
  44. F3.xin = [1,2,3]
  45. F3.KNOWN_BAD = {}
  46. def F4_powell(x):
  47. A = 1e4
  48. return [A*x[0]*x[1] - 1, np.exp(-x[0]) + np.exp(-x[1]) - (1 + 1/A)]
  49. F4_powell.xin = [-1, -2]
  50. F4_powell.KNOWN_BAD = {'linearmixing': nonlin.linearmixing,
  51. 'excitingmixing': nonlin.excitingmixing,
  52. 'diagbroyden': nonlin.diagbroyden}
  53. def F5(x):
  54. return pressure_network(x, 4, np.array([.5, .5, .5, .5]))
  55. F5.xin = [2., 0, 2, 0]
  56. F5.KNOWN_BAD = {'excitingmixing': nonlin.excitingmixing,
  57. 'linearmixing': nonlin.linearmixing,
  58. 'diagbroyden': nonlin.diagbroyden}
  59. def F6(x):
  60. x1, x2 = x
  61. J0 = np.array([[-4.256, 14.7],
  62. [0.8394989, 0.59964207]])
  63. v = np.array([(x1 + 3) * (x2**5 - 7) + 3*6,
  64. np.sin(x2 * np.exp(x1) - 1)])
  65. return -np.linalg.solve(J0, v)
  66. F6.xin = [-0.5, 1.4]
  67. F6.KNOWN_BAD = {'excitingmixing': nonlin.excitingmixing,
  68. 'linearmixing': nonlin.linearmixing,
  69. 'diagbroyden': nonlin.diagbroyden}
  70. #-------------------------------------------------------------------------------
  71. # Tests
  72. #-------------------------------------------------------------------------------
  73. class TestNonlin(object):
  74. """
  75. Check the Broyden methods for a few test problems.
  76. broyden1, broyden2, and newton_krylov must succeed for
  77. all functions. Some of the others don't -- tests in KNOWN_BAD are skipped.
  78. """
  79. def _check_nonlin_func(self, f, func, f_tol=1e-2):
  80. x = func(f, f.xin, f_tol=f_tol, maxiter=200, verbose=0)
  81. assert_(np.absolute(f(x)).max() < f_tol)
  82. def _check_root(self, f, method, f_tol=1e-2):
  83. res = root(f, f.xin, method=method,
  84. options={'ftol': f_tol, 'maxiter': 200, 'disp': 0})
  85. assert_(np.absolute(res.fun).max() < f_tol)
  86. @pytest.mark.xfail
  87. def _check_func_fail(self, *a, **kw):
  88. pass
  89. def test_problem_nonlin(self):
  90. for f in [F, F2, F2_lucky, F3, F4_powell, F5, F6]:
  91. for func in SOLVERS.values():
  92. if func in f.KNOWN_BAD.values():
  93. if func in MUST_WORK.values():
  94. self._check_func_fail(f, func)
  95. continue
  96. self._check_nonlin_func(f, func)
  97. def test_tol_norm_called(self):
  98. # Check that supplying tol_norm keyword to nonlin_solve works
  99. self._tol_norm_used = False
  100. def local_norm_func(x):
  101. self._tol_norm_used = True
  102. return np.absolute(x).max()
  103. nonlin.newton_krylov(F, F.xin, f_tol=1e-2, maxiter=200, verbose=0,
  104. tol_norm=local_norm_func)
  105. assert_(self._tol_norm_used)
  106. def test_problem_root(self):
  107. for f in [F, F2, F2_lucky, F3, F4_powell, F5, F6]:
  108. for meth in SOLVERS:
  109. if meth in f.KNOWN_BAD:
  110. if meth in MUST_WORK:
  111. self._check_func_fail(f, meth)
  112. continue
  113. self._check_root(f, meth)
  114. class TestSecant(object):
  115. """Check that some Jacobian approximations satisfy the secant condition"""
  116. xs = [np.array([1,2,3,4,5], float),
  117. np.array([2,3,4,5,1], float),
  118. np.array([3,4,5,1,2], float),
  119. np.array([4,5,1,2,3], float),
  120. np.array([9,1,9,1,3], float),
  121. np.array([0,1,9,1,3], float),
  122. np.array([5,5,7,1,1], float),
  123. np.array([1,2,7,5,1], float),]
  124. fs = [x**2 - 1 for x in xs]
  125. def _check_secant(self, jac_cls, npoints=1, **kw):
  126. """
  127. Check that the given Jacobian approximation satisfies secant
  128. conditions for last `npoints` points.
  129. """
  130. jac = jac_cls(**kw)
  131. jac.setup(self.xs[0], self.fs[0], None)
  132. for j, (x, f) in enumerate(zip(self.xs[1:], self.fs[1:])):
  133. jac.update(x, f)
  134. for k in xrange(min(npoints, j+1)):
  135. dx = self.xs[j-k+1] - self.xs[j-k]
  136. df = self.fs[j-k+1] - self.fs[j-k]
  137. assert_(np.allclose(dx, jac.solve(df)))
  138. # Check that the `npoints` secant bound is strict
  139. if j >= npoints:
  140. dx = self.xs[j-npoints+1] - self.xs[j-npoints]
  141. df = self.fs[j-npoints+1] - self.fs[j-npoints]
  142. assert_(not np.allclose(dx, jac.solve(df)))
  143. def test_broyden1(self):
  144. self._check_secant(nonlin.BroydenFirst)
  145. def test_broyden2(self):
  146. self._check_secant(nonlin.BroydenSecond)
  147. def test_broyden1_update(self):
  148. # Check that BroydenFirst update works as for a dense matrix
  149. jac = nonlin.BroydenFirst(alpha=0.1)
  150. jac.setup(self.xs[0], self.fs[0], None)
  151. B = np.identity(5) * (-1/0.1)
  152. for last_j, (x, f) in enumerate(zip(self.xs[1:], self.fs[1:])):
  153. df = f - self.fs[last_j]
  154. dx = x - self.xs[last_j]
  155. B += (df - dot(B, dx))[:,None] * dx[None,:] / dot(dx, dx)
  156. jac.update(x, f)
  157. assert_(np.allclose(jac.todense(), B, rtol=1e-10, atol=1e-13))
  158. def test_broyden2_update(self):
  159. # Check that BroydenSecond update works as for a dense matrix
  160. jac = nonlin.BroydenSecond(alpha=0.1)
  161. jac.setup(self.xs[0], self.fs[0], None)
  162. H = np.identity(5) * (-0.1)
  163. for last_j, (x, f) in enumerate(zip(self.xs[1:], self.fs[1:])):
  164. df = f - self.fs[last_j]
  165. dx = x - self.xs[last_j]
  166. H += (dx - dot(H, df))[:,None] * df[None,:] / dot(df, df)
  167. jac.update(x, f)
  168. assert_(np.allclose(jac.todense(), inv(H), rtol=1e-10, atol=1e-13))
  169. def test_anderson(self):
  170. # Anderson mixing (with w0=0) satisfies secant conditions
  171. # for the last M iterates, see [Ey]_
  172. #
  173. # .. [Ey] V. Eyert, J. Comp. Phys., 124, 271 (1996).
  174. self._check_secant(nonlin.Anderson, M=3, w0=0, npoints=3)
  175. class TestLinear(object):
  176. """Solve a linear equation;
  177. some methods find the exact solution in a finite number of steps"""
  178. def _check(self, jac, N, maxiter, complex=False, **kw):
  179. np.random.seed(123)
  180. A = np.random.randn(N, N)
  181. if complex:
  182. A = A + 1j*np.random.randn(N, N)
  183. b = np.random.randn(N)
  184. if complex:
  185. b = b + 1j*np.random.randn(N)
  186. def func(x):
  187. return dot(A, x) - b
  188. sol = nonlin.nonlin_solve(func, np.zeros(N), jac, maxiter=maxiter,
  189. f_tol=1e-6, line_search=None, verbose=0)
  190. assert_(np.allclose(dot(A, sol), b, atol=1e-6))
  191. def test_broyden1(self):
  192. # Broyden methods solve linear systems exactly in 2*N steps
  193. self._check(nonlin.BroydenFirst(alpha=1.0), 20, 41, False)
  194. self._check(nonlin.BroydenFirst(alpha=1.0), 20, 41, True)
  195. def test_broyden2(self):
  196. # Broyden methods solve linear systems exactly in 2*N steps
  197. self._check(nonlin.BroydenSecond(alpha=1.0), 20, 41, False)
  198. self._check(nonlin.BroydenSecond(alpha=1.0), 20, 41, True)
  199. def test_anderson(self):
  200. # Anderson is rather similar to Broyden, if given enough storage space
  201. self._check(nonlin.Anderson(M=50, alpha=1.0), 20, 29, False)
  202. self._check(nonlin.Anderson(M=50, alpha=1.0), 20, 29, True)
  203. def test_krylov(self):
  204. # Krylov methods solve linear systems exactly in N inner steps
  205. self._check(nonlin.KrylovJacobian, 20, 2, False, inner_m=10)
  206. self._check(nonlin.KrylovJacobian, 20, 2, True, inner_m=10)
  207. class TestJacobianDotSolve(object):
  208. """Check that solve/dot methods in Jacobian approximations are consistent"""
  209. def _func(self, x):
  210. return x**2 - 1 + np.dot(self.A, x)
  211. def _check_dot(self, jac_cls, complex=False, tol=1e-6, **kw):
  212. np.random.seed(123)
  213. N = 7
  214. def rand(*a):
  215. q = np.random.rand(*a)
  216. if complex:
  217. q = q + 1j*np.random.rand(*a)
  218. return q
  219. def assert_close(a, b, msg):
  220. d = abs(a - b).max()
  221. f = tol + abs(b).max()*tol
  222. if d > f:
  223. raise AssertionError('%s: err %g' % (msg, d))
  224. self.A = rand(N, N)
  225. # initialize
  226. x0 = np.random.rand(N)
  227. jac = jac_cls(**kw)
  228. jac.setup(x0, self._func(x0), self._func)
  229. # check consistency
  230. for k in xrange(2*N):
  231. v = rand(N)
  232. if hasattr(jac, '__array__'):
  233. Jd = np.array(jac)
  234. if hasattr(jac, 'solve'):
  235. Gv = jac.solve(v)
  236. Gv2 = np.linalg.solve(Jd, v)
  237. assert_close(Gv, Gv2, 'solve vs array')
  238. if hasattr(jac, 'rsolve'):
  239. Gv = jac.rsolve(v)
  240. Gv2 = np.linalg.solve(Jd.T.conj(), v)
  241. assert_close(Gv, Gv2, 'rsolve vs array')
  242. if hasattr(jac, 'matvec'):
  243. Jv = jac.matvec(v)
  244. Jv2 = np.dot(Jd, v)
  245. assert_close(Jv, Jv2, 'dot vs array')
  246. if hasattr(jac, 'rmatvec'):
  247. Jv = jac.rmatvec(v)
  248. Jv2 = np.dot(Jd.T.conj(), v)
  249. assert_close(Jv, Jv2, 'rmatvec vs array')
  250. if hasattr(jac, 'matvec') and hasattr(jac, 'solve'):
  251. Jv = jac.matvec(v)
  252. Jv2 = jac.solve(jac.matvec(Jv))
  253. assert_close(Jv, Jv2, 'dot vs solve')
  254. if hasattr(jac, 'rmatvec') and hasattr(jac, 'rsolve'):
  255. Jv = jac.rmatvec(v)
  256. Jv2 = jac.rmatvec(jac.rsolve(Jv))
  257. assert_close(Jv, Jv2, 'rmatvec vs rsolve')
  258. x = rand(N)
  259. jac.update(x, self._func(x))
  260. def test_broyden1(self):
  261. self._check_dot(nonlin.BroydenFirst, complex=False)
  262. self._check_dot(nonlin.BroydenFirst, complex=True)
  263. def test_broyden2(self):
  264. self._check_dot(nonlin.BroydenSecond, complex=False)
  265. self._check_dot(nonlin.BroydenSecond, complex=True)
  266. def test_anderson(self):
  267. self._check_dot(nonlin.Anderson, complex=False)
  268. self._check_dot(nonlin.Anderson, complex=True)
  269. def test_diagbroyden(self):
  270. self._check_dot(nonlin.DiagBroyden, complex=False)
  271. self._check_dot(nonlin.DiagBroyden, complex=True)
  272. def test_linearmixing(self):
  273. self._check_dot(nonlin.LinearMixing, complex=False)
  274. self._check_dot(nonlin.LinearMixing, complex=True)
  275. def test_excitingmixing(self):
  276. self._check_dot(nonlin.ExcitingMixing, complex=False)
  277. self._check_dot(nonlin.ExcitingMixing, complex=True)
  278. def test_krylov(self):
  279. self._check_dot(nonlin.KrylovJacobian, complex=False, tol=1e-3)
  280. self._check_dot(nonlin.KrylovJacobian, complex=True, tol=1e-3)
  281. class TestNonlinOldTests(object):
  282. """ Test case for a simple constrained entropy maximization problem
  283. (the machine translation example of Berger et al in
  284. Computational Linguistics, vol 22, num 1, pp 39--72, 1996.)
  285. """
  286. def test_broyden1(self):
  287. x = nonlin.broyden1(F,F.xin,iter=12,alpha=1)
  288. assert_(nonlin.norm(x) < 1e-9)
  289. assert_(nonlin.norm(F(x)) < 1e-9)
  290. def test_broyden2(self):
  291. x = nonlin.broyden2(F,F.xin,iter=12,alpha=1)
  292. assert_(nonlin.norm(x) < 1e-9)
  293. assert_(nonlin.norm(F(x)) < 1e-9)
  294. def test_anderson(self):
  295. x = nonlin.anderson(F,F.xin,iter=12,alpha=0.03,M=5)
  296. assert_(nonlin.norm(x) < 0.33)
  297. def test_linearmixing(self):
  298. x = nonlin.linearmixing(F,F.xin,iter=60,alpha=0.5)
  299. assert_(nonlin.norm(x) < 1e-7)
  300. assert_(nonlin.norm(F(x)) < 1e-7)
  301. def test_exciting(self):
  302. x = nonlin.excitingmixing(F,F.xin,iter=20,alpha=0.5)
  303. assert_(nonlin.norm(x) < 1e-5)
  304. assert_(nonlin.norm(F(x)) < 1e-5)
  305. def test_diagbroyden(self):
  306. x = nonlin.diagbroyden(F,F.xin,iter=11,alpha=1)
  307. assert_(nonlin.norm(x) < 1e-8)
  308. assert_(nonlin.norm(F(x)) < 1e-8)
  309. def test_root_broyden1(self):
  310. res = root(F, F.xin, method='broyden1',
  311. options={'nit': 12, 'jac_options': {'alpha': 1}})
  312. assert_(nonlin.norm(res.x) < 1e-9)
  313. assert_(nonlin.norm(res.fun) < 1e-9)
  314. def test_root_broyden2(self):
  315. res = root(F, F.xin, method='broyden2',
  316. options={'nit': 12, 'jac_options': {'alpha': 1}})
  317. assert_(nonlin.norm(res.x) < 1e-9)
  318. assert_(nonlin.norm(res.fun) < 1e-9)
  319. def test_root_anderson(self):
  320. res = root(F, F.xin, method='anderson',
  321. options={'nit': 12,
  322. 'jac_options': {'alpha': 0.03, 'M': 5}})
  323. assert_(nonlin.norm(res.x) < 0.33)
  324. def test_root_linearmixing(self):
  325. res = root(F, F.xin, method='linearmixing',
  326. options={'nit': 60,
  327. 'jac_options': {'alpha': 0.5}})
  328. assert_(nonlin.norm(res.x) < 1e-7)
  329. assert_(nonlin.norm(res.fun) < 1e-7)
  330. def test_root_excitingmixing(self):
  331. res = root(F, F.xin, method='excitingmixing',
  332. options={'nit': 20,
  333. 'jac_options': {'alpha': 0.5}})
  334. assert_(nonlin.norm(res.x) < 1e-5)
  335. assert_(nonlin.norm(res.fun) < 1e-5)
  336. def test_root_diagbroyden(self):
  337. res = root(F, F.xin, method='diagbroyden',
  338. options={'nit': 11,
  339. 'jac_options': {'alpha': 1}})
  340. assert_(nonlin.norm(res.x) < 1e-8)
  341. assert_(nonlin.norm(res.fun) < 1e-8)