test_differentiable_functions.py 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731
  1. import pytest
  2. import numpy as np
  3. from numpy.testing import (TestCase, assert_array_almost_equal,
  4. assert_array_equal, assert_, assert_allclose,
  5. assert_equal)
  6. from scipy.sparse import csr_matrix
  7. from scipy.sparse.linalg import LinearOperator
  8. from scipy.optimize._differentiable_functions import (ScalarFunction,
  9. VectorFunction,
  10. LinearVectorFunction,
  11. IdentityVectorFunction)
  12. from scipy.optimize import rosen, rosen_der, rosen_hess
  13. from scipy.optimize._hessian_update_strategy import BFGS
  14. class ExScalarFunction:
  15. def __init__(self):
  16. self.nfev = 0
  17. self.ngev = 0
  18. self.nhev = 0
  19. def fun(self, x):
  20. self.nfev += 1
  21. return 2*(x[0]**2 + x[1]**2 - 1) - x[0]
  22. def grad(self, x):
  23. self.ngev += 1
  24. return np.array([4*x[0]-1, 4*x[1]])
  25. def hess(self, x):
  26. self.nhev += 1
  27. return 4*np.eye(2)
  28. class TestScalarFunction(TestCase):
  29. def test_finite_difference_grad(self):
  30. ex = ExScalarFunction()
  31. nfev = 0
  32. ngev = 0
  33. x0 = [1.0, 0.0]
  34. analit = ScalarFunction(ex.fun, x0, (), ex.grad,
  35. ex.hess, None, (-np.inf, np.inf))
  36. nfev += 1
  37. ngev += 1
  38. assert_array_equal(ex.nfev, nfev)
  39. assert_array_equal(analit.nfev, nfev)
  40. assert_array_equal(ex.ngev, ngev)
  41. assert_array_equal(analit.ngev, nfev)
  42. approx = ScalarFunction(ex.fun, x0, (), '2-point',
  43. ex.hess, None, (-np.inf, np.inf))
  44. nfev += 3
  45. ngev += 1
  46. assert_array_equal(ex.nfev, nfev)
  47. assert_array_equal(analit.nfev+approx.nfev, nfev)
  48. assert_array_equal(analit.ngev+approx.ngev, ngev)
  49. assert_array_equal(analit.f, approx.f)
  50. assert_array_almost_equal(analit.g, approx.g)
  51. x = [10, 0.3]
  52. f_analit = analit.fun(x)
  53. g_analit = analit.grad(x)
  54. nfev += 1
  55. ngev += 1
  56. assert_array_equal(ex.nfev, nfev)
  57. assert_array_equal(analit.nfev+approx.nfev, nfev)
  58. assert_array_equal(analit.ngev+approx.ngev, ngev)
  59. f_approx = approx.fun(x)
  60. g_approx = approx.grad(x)
  61. nfev += 3
  62. ngev += 1
  63. assert_array_equal(ex.nfev, nfev)
  64. assert_array_equal(analit.nfev+approx.nfev, nfev)
  65. assert_array_equal(analit.ngev+approx.ngev, ngev)
  66. assert_array_almost_equal(f_analit, f_approx)
  67. assert_array_almost_equal(g_analit, g_approx)
  68. x = [2.0, 1.0]
  69. g_analit = analit.grad(x)
  70. ngev += 1
  71. assert_array_equal(ex.nfev, nfev)
  72. assert_array_equal(analit.nfev+approx.nfev, nfev)
  73. assert_array_equal(analit.ngev+approx.ngev, ngev)
  74. g_approx = approx.grad(x)
  75. nfev += 3
  76. ngev += 1
  77. assert_array_equal(ex.nfev, nfev)
  78. assert_array_equal(analit.nfev+approx.nfev, nfev)
  79. assert_array_equal(analit.ngev+approx.ngev, ngev)
  80. assert_array_almost_equal(g_analit, g_approx)
  81. x = [2.5, 0.3]
  82. f_analit = analit.fun(x)
  83. g_analit = analit.grad(x)
  84. nfev += 1
  85. ngev += 1
  86. assert_array_equal(ex.nfev, nfev)
  87. assert_array_equal(analit.nfev+approx.nfev, nfev)
  88. assert_array_equal(analit.ngev+approx.ngev, ngev)
  89. f_approx = approx.fun(x)
  90. g_approx = approx.grad(x)
  91. nfev += 3
  92. ngev += 1
  93. assert_array_equal(ex.nfev, nfev)
  94. assert_array_equal(analit.nfev+approx.nfev, nfev)
  95. assert_array_equal(analit.ngev+approx.ngev, ngev)
  96. assert_array_almost_equal(f_analit, f_approx)
  97. assert_array_almost_equal(g_analit, g_approx)
  98. x = [2, 0.3]
  99. f_analit = analit.fun(x)
  100. g_analit = analit.grad(x)
  101. nfev += 1
  102. ngev += 1
  103. assert_array_equal(ex.nfev, nfev)
  104. assert_array_equal(analit.nfev+approx.nfev, nfev)
  105. assert_array_equal(analit.ngev+approx.ngev, ngev)
  106. f_approx = approx.fun(x)
  107. g_approx = approx.grad(x)
  108. nfev += 3
  109. ngev += 1
  110. assert_array_equal(ex.nfev, nfev)
  111. assert_array_equal(analit.nfev+approx.nfev, nfev)
  112. assert_array_equal(analit.ngev+approx.ngev, ngev)
  113. assert_array_almost_equal(f_analit, f_approx)
  114. assert_array_almost_equal(g_analit, g_approx)
  115. def test_fun_and_grad(self):
  116. ex = ExScalarFunction()
  117. def fg_allclose(x, y):
  118. assert_allclose(x[0], y[0])
  119. assert_allclose(x[1], y[1])
  120. # with analytic gradient
  121. x0 = [2.0, 0.3]
  122. analit = ScalarFunction(ex.fun, x0, (), ex.grad,
  123. ex.hess, None, (-np.inf, np.inf))
  124. fg = ex.fun(x0), ex.grad(x0)
  125. fg_allclose(analit.fun_and_grad(x0), fg)
  126. assert analit.ngev == 1
  127. x0[1] = 1.
  128. fg = ex.fun(x0), ex.grad(x0)
  129. fg_allclose(analit.fun_and_grad(x0), fg)
  130. # with finite difference gradient
  131. x0 = [2.0, 0.3]
  132. sf = ScalarFunction(ex.fun, x0, (), '3-point',
  133. ex.hess, None, (-np.inf, np.inf))
  134. assert sf.ngev == 1
  135. fg = ex.fun(x0), ex.grad(x0)
  136. fg_allclose(sf.fun_and_grad(x0), fg)
  137. assert sf.ngev == 1
  138. x0[1] = 1.
  139. fg = ex.fun(x0), ex.grad(x0)
  140. fg_allclose(sf.fun_and_grad(x0), fg)
  141. def test_finite_difference_hess_linear_operator(self):
  142. ex = ExScalarFunction()
  143. nfev = 0
  144. ngev = 0
  145. nhev = 0
  146. x0 = [1.0, 0.0]
  147. analit = ScalarFunction(ex.fun, x0, (), ex.grad,
  148. ex.hess, None, (-np.inf, np.inf))
  149. nfev += 1
  150. ngev += 1
  151. nhev += 1
  152. assert_array_equal(ex.nfev, nfev)
  153. assert_array_equal(analit.nfev, nfev)
  154. assert_array_equal(ex.ngev, ngev)
  155. assert_array_equal(analit.ngev, ngev)
  156. assert_array_equal(ex.nhev, nhev)
  157. assert_array_equal(analit.nhev, nhev)
  158. approx = ScalarFunction(ex.fun, x0, (), ex.grad,
  159. '2-point', None, (-np.inf, np.inf))
  160. assert_(isinstance(approx.H, LinearOperator))
  161. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  162. assert_array_equal(analit.f, approx.f)
  163. assert_array_almost_equal(analit.g, approx.g)
  164. assert_array_almost_equal(analit.H.dot(v), approx.H.dot(v))
  165. nfev += 1
  166. ngev += 4
  167. assert_array_equal(ex.nfev, nfev)
  168. assert_array_equal(analit.nfev+approx.nfev, nfev)
  169. assert_array_equal(ex.ngev, ngev)
  170. assert_array_equal(analit.ngev+approx.ngev, ngev)
  171. assert_array_equal(ex.nhev, nhev)
  172. assert_array_equal(analit.nhev+approx.nhev, nhev)
  173. x = [2.0, 1.0]
  174. H_analit = analit.hess(x)
  175. nhev += 1
  176. assert_array_equal(ex.nfev, nfev)
  177. assert_array_equal(analit.nfev+approx.nfev, nfev)
  178. assert_array_equal(ex.ngev, ngev)
  179. assert_array_equal(analit.ngev+approx.ngev, ngev)
  180. assert_array_equal(ex.nhev, nhev)
  181. assert_array_equal(analit.nhev+approx.nhev, nhev)
  182. H_approx = approx.hess(x)
  183. assert_(isinstance(H_approx, LinearOperator))
  184. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  185. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v))
  186. ngev += 4
  187. assert_array_equal(ex.nfev, nfev)
  188. assert_array_equal(analit.nfev+approx.nfev, nfev)
  189. assert_array_equal(ex.ngev, ngev)
  190. assert_array_equal(analit.ngev+approx.ngev, ngev)
  191. assert_array_equal(ex.nhev, nhev)
  192. assert_array_equal(analit.nhev+approx.nhev, nhev)
  193. x = [2.1, 1.2]
  194. H_analit = analit.hess(x)
  195. nhev += 1
  196. assert_array_equal(ex.nfev, nfev)
  197. assert_array_equal(analit.nfev+approx.nfev, nfev)
  198. assert_array_equal(ex.ngev, ngev)
  199. assert_array_equal(analit.ngev+approx.ngev, ngev)
  200. assert_array_equal(ex.nhev, nhev)
  201. assert_array_equal(analit.nhev+approx.nhev, nhev)
  202. H_approx = approx.hess(x)
  203. assert_(isinstance(H_approx, LinearOperator))
  204. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  205. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v))
  206. ngev += 4
  207. assert_array_equal(ex.nfev, nfev)
  208. assert_array_equal(analit.nfev+approx.nfev, nfev)
  209. assert_array_equal(ex.ngev, ngev)
  210. assert_array_equal(analit.ngev+approx.ngev, ngev)
  211. assert_array_equal(ex.nhev, nhev)
  212. assert_array_equal(analit.nhev+approx.nhev, nhev)
  213. x = [2.5, 0.3]
  214. _ = analit.grad(x)
  215. H_analit = analit.hess(x)
  216. ngev += 1
  217. nhev += 1
  218. assert_array_equal(ex.nfev, nfev)
  219. assert_array_equal(analit.nfev+approx.nfev, nfev)
  220. assert_array_equal(ex.ngev, ngev)
  221. assert_array_equal(analit.ngev+approx.ngev, ngev)
  222. assert_array_equal(ex.nhev, nhev)
  223. assert_array_equal(analit.nhev+approx.nhev, nhev)
  224. _ = approx.grad(x)
  225. H_approx = approx.hess(x)
  226. assert_(isinstance(H_approx, LinearOperator))
  227. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  228. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v))
  229. ngev += 4
  230. assert_array_equal(ex.nfev, nfev)
  231. assert_array_equal(analit.nfev+approx.nfev, nfev)
  232. assert_array_equal(ex.ngev, ngev)
  233. assert_array_equal(analit.ngev+approx.ngev, ngev)
  234. assert_array_equal(ex.nhev, nhev)
  235. assert_array_equal(analit.nhev+approx.nhev, nhev)
  236. x = [5.2, 2.3]
  237. _ = analit.grad(x)
  238. H_analit = analit.hess(x)
  239. ngev += 1
  240. nhev += 1
  241. assert_array_equal(ex.nfev, nfev)
  242. assert_array_equal(analit.nfev+approx.nfev, nfev)
  243. assert_array_equal(ex.ngev, ngev)
  244. assert_array_equal(analit.ngev+approx.ngev, ngev)
  245. assert_array_equal(ex.nhev, nhev)
  246. assert_array_equal(analit.nhev+approx.nhev, nhev)
  247. _ = approx.grad(x)
  248. H_approx = approx.hess(x)
  249. assert_(isinstance(H_approx, LinearOperator))
  250. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  251. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v))
  252. ngev += 4
  253. assert_array_equal(ex.nfev, nfev)
  254. assert_array_equal(analit.nfev+approx.nfev, nfev)
  255. assert_array_equal(ex.ngev, ngev)
  256. assert_array_equal(analit.ngev+approx.ngev, ngev)
  257. assert_array_equal(ex.nhev, nhev)
  258. assert_array_equal(analit.nhev+approx.nhev, nhev)
  259. def test_x_storage_overlap(self):
  260. # Scalar_Function should not store references to arrays, it should
  261. # store copies - this checks that updating an array in-place causes
  262. # Scalar_Function.x to be updated.
  263. def f(x):
  264. return np.sum(np.asarray(x) ** 2)
  265. x = np.array([1., 2., 3.])
  266. sf = ScalarFunction(f, x, (), '3-point', lambda x: x, None, (-np.inf, np.inf))
  267. assert x is not sf.x
  268. assert_equal(sf.fun(x), 14.0)
  269. assert x is not sf.x
  270. x[0] = 0.
  271. f1 = sf.fun(x)
  272. assert_equal(f1, 13.0)
  273. x[0] = 1
  274. f2 = sf.fun(x)
  275. assert_equal(f2, 14.0)
  276. assert x is not sf.x
  277. # now test with a HessianUpdate strategy specified
  278. hess = BFGS()
  279. x = np.array([1., 2., 3.])
  280. sf = ScalarFunction(f, x, (), '3-point', hess, None, (-np.inf, np.inf))
  281. assert x is not sf.x
  282. assert_equal(sf.fun(x), 14.0)
  283. assert x is not sf.x
  284. x[0] = 0.
  285. f1 = sf.fun(x)
  286. assert_equal(f1, 13.0)
  287. x[0] = 1
  288. f2 = sf.fun(x)
  289. assert_equal(f2, 14.0)
  290. assert x is not sf.x
  291. # gh13740 x is changed in user function
  292. def ff(x):
  293. x *= x # overwrite x
  294. return np.sum(x)
  295. x = np.array([1., 2., 3.])
  296. sf = ScalarFunction(
  297. ff, x, (), '3-point', lambda x: x, None, (-np.inf, np.inf)
  298. )
  299. assert x is not sf.x
  300. assert_equal(sf.fun(x), 14.0)
  301. assert_equal(sf.x, np.array([1., 2., 3.]))
  302. assert x is not sf.x
  303. def test_lowest_x(self):
  304. # ScalarFunction should remember the lowest func(x) visited.
  305. x0 = np.array([2, 3, 4])
  306. sf = ScalarFunction(rosen, x0, (), rosen_der, rosen_hess,
  307. None, None)
  308. sf.fun([1, 1, 1])
  309. sf.fun(x0)
  310. sf.fun([1.01, 1, 1.0])
  311. sf.grad([1.01, 1, 1.0])
  312. assert_equal(sf._lowest_f, 0.0)
  313. assert_equal(sf._lowest_x, [1.0, 1.0, 1.0])
  314. sf = ScalarFunction(rosen, x0, (), '2-point', rosen_hess,
  315. None, (-np.inf, np.inf))
  316. sf.fun([1, 1, 1])
  317. sf.fun(x0)
  318. sf.fun([1.01, 1, 1.0])
  319. sf.grad([1.01, 1, 1.0])
  320. assert_equal(sf._lowest_f, 0.0)
  321. assert_equal(sf._lowest_x, [1.0, 1.0, 1.0])
  322. class ExVectorialFunction:
  323. def __init__(self):
  324. self.nfev = 0
  325. self.njev = 0
  326. self.nhev = 0
  327. def fun(self, x):
  328. self.nfev += 1
  329. return np.array([2*(x[0]**2 + x[1]**2 - 1) - x[0],
  330. 4*(x[0]**3 + x[1]**2 - 4) - 3*x[0]])
  331. def jac(self, x):
  332. self.njev += 1
  333. return np.array([[4*x[0]-1, 4*x[1]],
  334. [12*x[0]**2-3, 8*x[1]]])
  335. def hess(self, x, v):
  336. self.nhev += 1
  337. return v[0]*4*np.eye(2) + v[1]*np.array([[24*x[0], 0],
  338. [0, 8]])
  339. class TestVectorialFunction(TestCase):
  340. def test_finite_difference_jac(self):
  341. ex = ExVectorialFunction()
  342. nfev = 0
  343. njev = 0
  344. x0 = [1.0, 0.0]
  345. analit = VectorFunction(ex.fun, x0, ex.jac, ex.hess, None, None,
  346. (-np.inf, np.inf), None)
  347. nfev += 1
  348. njev += 1
  349. assert_array_equal(ex.nfev, nfev)
  350. assert_array_equal(analit.nfev, nfev)
  351. assert_array_equal(ex.njev, njev)
  352. assert_array_equal(analit.njev, njev)
  353. approx = VectorFunction(ex.fun, x0, '2-point', ex.hess, None, None,
  354. (-np.inf, np.inf), None)
  355. nfev += 3
  356. assert_array_equal(ex.nfev, nfev)
  357. assert_array_equal(analit.nfev+approx.nfev, nfev)
  358. assert_array_equal(ex.njev, njev)
  359. assert_array_equal(analit.njev+approx.njev, njev)
  360. assert_array_equal(analit.f, approx.f)
  361. assert_array_almost_equal(analit.J, approx.J)
  362. x = [10, 0.3]
  363. f_analit = analit.fun(x)
  364. J_analit = analit.jac(x)
  365. nfev += 1
  366. njev += 1
  367. assert_array_equal(ex.nfev, nfev)
  368. assert_array_equal(analit.nfev+approx.nfev, nfev)
  369. assert_array_equal(ex.njev, njev)
  370. assert_array_equal(analit.njev+approx.njev, njev)
  371. f_approx = approx.fun(x)
  372. J_approx = approx.jac(x)
  373. nfev += 3
  374. assert_array_equal(ex.nfev, nfev)
  375. assert_array_equal(analit.nfev+approx.nfev, nfev)
  376. assert_array_equal(ex.njev, njev)
  377. assert_array_equal(analit.njev+approx.njev, njev)
  378. assert_array_almost_equal(f_analit, f_approx)
  379. assert_array_almost_equal(J_analit, J_approx, decimal=4)
  380. x = [2.0, 1.0]
  381. J_analit = analit.jac(x)
  382. njev += 1
  383. assert_array_equal(ex.nfev, nfev)
  384. assert_array_equal(analit.nfev+approx.nfev, nfev)
  385. assert_array_equal(ex.njev, njev)
  386. assert_array_equal(analit.njev+approx.njev, njev)
  387. J_approx = approx.jac(x)
  388. nfev += 3
  389. assert_array_equal(ex.nfev, nfev)
  390. assert_array_equal(analit.nfev+approx.nfev, nfev)
  391. assert_array_equal(ex.njev, njev)
  392. assert_array_equal(analit.njev+approx.njev, njev)
  393. assert_array_almost_equal(J_analit, J_approx)
  394. x = [2.5, 0.3]
  395. f_analit = analit.fun(x)
  396. J_analit = analit.jac(x)
  397. nfev += 1
  398. njev += 1
  399. assert_array_equal(ex.nfev, nfev)
  400. assert_array_equal(analit.nfev+approx.nfev, nfev)
  401. assert_array_equal(ex.njev, njev)
  402. assert_array_equal(analit.njev+approx.njev, njev)
  403. f_approx = approx.fun(x)
  404. J_approx = approx.jac(x)
  405. nfev += 3
  406. assert_array_equal(ex.nfev, nfev)
  407. assert_array_equal(analit.nfev+approx.nfev, nfev)
  408. assert_array_equal(ex.njev, njev)
  409. assert_array_equal(analit.njev+approx.njev, njev)
  410. assert_array_almost_equal(f_analit, f_approx)
  411. assert_array_almost_equal(J_analit, J_approx)
  412. x = [2, 0.3]
  413. f_analit = analit.fun(x)
  414. J_analit = analit.jac(x)
  415. nfev += 1
  416. njev += 1
  417. assert_array_equal(ex.nfev, nfev)
  418. assert_array_equal(analit.nfev+approx.nfev, nfev)
  419. assert_array_equal(ex.njev, njev)
  420. assert_array_equal(analit.njev+approx.njev, njev)
  421. f_approx = approx.fun(x)
  422. J_approx = approx.jac(x)
  423. nfev += 3
  424. assert_array_equal(ex.nfev, nfev)
  425. assert_array_equal(analit.nfev+approx.nfev, nfev)
  426. assert_array_equal(ex.njev, njev)
  427. assert_array_equal(analit.njev+approx.njev, njev)
  428. assert_array_almost_equal(f_analit, f_approx)
  429. assert_array_almost_equal(J_analit, J_approx)
  430. def test_finite_difference_hess_linear_operator(self):
  431. ex = ExVectorialFunction()
  432. nfev = 0
  433. njev = 0
  434. nhev = 0
  435. x0 = [1.0, 0.0]
  436. v0 = [1.0, 2.0]
  437. analit = VectorFunction(ex.fun, x0, ex.jac, ex.hess, None, None,
  438. (-np.inf, np.inf), None)
  439. nfev += 1
  440. njev += 1
  441. nhev += 1
  442. assert_array_equal(ex.nfev, nfev)
  443. assert_array_equal(analit.nfev, nfev)
  444. assert_array_equal(ex.njev, njev)
  445. assert_array_equal(analit.njev, njev)
  446. assert_array_equal(ex.nhev, nhev)
  447. assert_array_equal(analit.nhev, nhev)
  448. approx = VectorFunction(ex.fun, x0, ex.jac, '2-point', None, None,
  449. (-np.inf, np.inf), None)
  450. assert_(isinstance(approx.H, LinearOperator))
  451. for p in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  452. assert_array_equal(analit.f, approx.f)
  453. assert_array_almost_equal(analit.J, approx.J)
  454. assert_array_almost_equal(analit.H.dot(p), approx.H.dot(p))
  455. nfev += 1
  456. njev += 4
  457. assert_array_equal(ex.nfev, nfev)
  458. assert_array_equal(analit.nfev+approx.nfev, nfev)
  459. assert_array_equal(ex.njev, njev)
  460. assert_array_equal(analit.njev+approx.njev, njev)
  461. assert_array_equal(ex.nhev, nhev)
  462. assert_array_equal(analit.nhev+approx.nhev, nhev)
  463. x = [2.0, 1.0]
  464. H_analit = analit.hess(x, v0)
  465. nhev += 1
  466. assert_array_equal(ex.nfev, nfev)
  467. assert_array_equal(analit.nfev+approx.nfev, nfev)
  468. assert_array_equal(ex.njev, njev)
  469. assert_array_equal(analit.njev+approx.njev, njev)
  470. assert_array_equal(ex.nhev, nhev)
  471. assert_array_equal(analit.nhev+approx.nhev, nhev)
  472. H_approx = approx.hess(x, v0)
  473. assert_(isinstance(H_approx, LinearOperator))
  474. for p in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  475. assert_array_almost_equal(H_analit.dot(p), H_approx.dot(p),
  476. decimal=5)
  477. njev += 4
  478. assert_array_equal(ex.nfev, nfev)
  479. assert_array_equal(analit.nfev+approx.nfev, nfev)
  480. assert_array_equal(ex.njev, njev)
  481. assert_array_equal(analit.njev+approx.njev, njev)
  482. assert_array_equal(ex.nhev, nhev)
  483. assert_array_equal(analit.nhev+approx.nhev, nhev)
  484. x = [2.1, 1.2]
  485. v = [1.0, 1.0]
  486. H_analit = analit.hess(x, v)
  487. nhev += 1
  488. assert_array_equal(ex.nfev, nfev)
  489. assert_array_equal(analit.nfev+approx.nfev, nfev)
  490. assert_array_equal(ex.njev, njev)
  491. assert_array_equal(analit.njev+approx.njev, njev)
  492. assert_array_equal(ex.nhev, nhev)
  493. assert_array_equal(analit.nhev+approx.nhev, nhev)
  494. H_approx = approx.hess(x, v)
  495. assert_(isinstance(H_approx, LinearOperator))
  496. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  497. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v))
  498. njev += 4
  499. assert_array_equal(ex.nfev, nfev)
  500. assert_array_equal(analit.nfev+approx.nfev, nfev)
  501. assert_array_equal(ex.njev, njev)
  502. assert_array_equal(analit.njev+approx.njev, njev)
  503. assert_array_equal(ex.nhev, nhev)
  504. assert_array_equal(analit.nhev+approx.nhev, nhev)
  505. x = [2.5, 0.3]
  506. _ = analit.jac(x)
  507. H_analit = analit.hess(x, v0)
  508. njev += 1
  509. nhev += 1
  510. assert_array_equal(ex.nfev, nfev)
  511. assert_array_equal(analit.nfev+approx.nfev, nfev)
  512. assert_array_equal(ex.njev, njev)
  513. assert_array_equal(analit.njev+approx.njev, njev)
  514. assert_array_equal(ex.nhev, nhev)
  515. assert_array_equal(analit.nhev+approx.nhev, nhev)
  516. _ = approx.jac(x)
  517. H_approx = approx.hess(x, v0)
  518. assert_(isinstance(H_approx, LinearOperator))
  519. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  520. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v), decimal=4)
  521. njev += 4
  522. assert_array_equal(ex.nfev, nfev)
  523. assert_array_equal(analit.nfev+approx.nfev, nfev)
  524. assert_array_equal(ex.njev, njev)
  525. assert_array_equal(analit.njev+approx.njev, njev)
  526. assert_array_equal(ex.nhev, nhev)
  527. assert_array_equal(analit.nhev+approx.nhev, nhev)
  528. x = [5.2, 2.3]
  529. v = [2.3, 5.2]
  530. _ = analit.jac(x)
  531. H_analit = analit.hess(x, v)
  532. njev += 1
  533. nhev += 1
  534. assert_array_equal(ex.nfev, nfev)
  535. assert_array_equal(analit.nfev+approx.nfev, nfev)
  536. assert_array_equal(ex.njev, njev)
  537. assert_array_equal(analit.njev+approx.njev, njev)
  538. assert_array_equal(ex.nhev, nhev)
  539. assert_array_equal(analit.nhev+approx.nhev, nhev)
  540. _ = approx.jac(x)
  541. H_approx = approx.hess(x, v)
  542. assert_(isinstance(H_approx, LinearOperator))
  543. for v in ([1.0, 2.0], [3.0, 4.0], [5.0, 2.0]):
  544. assert_array_almost_equal(H_analit.dot(v), H_approx.dot(v), decimal=4)
  545. njev += 4
  546. assert_array_equal(ex.nfev, nfev)
  547. assert_array_equal(analit.nfev+approx.nfev, nfev)
  548. assert_array_equal(ex.njev, njev)
  549. assert_array_equal(analit.njev+approx.njev, njev)
  550. assert_array_equal(ex.nhev, nhev)
  551. assert_array_equal(analit.nhev+approx.nhev, nhev)
  552. def test_x_storage_overlap(self):
  553. # VectorFunction should not store references to arrays, it should
  554. # store copies - this checks that updating an array in-place causes
  555. # Scalar_Function.x to be updated.
  556. ex = ExVectorialFunction()
  557. x0 = np.array([1.0, 0.0])
  558. vf = VectorFunction(ex.fun, x0, '3-point', ex.hess, None, None,
  559. (-np.inf, np.inf), None)
  560. assert x0 is not vf.x
  561. assert_equal(vf.fun(x0), ex.fun(x0))
  562. assert x0 is not vf.x
  563. x0[0] = 2.
  564. assert_equal(vf.fun(x0), ex.fun(x0))
  565. assert x0 is not vf.x
  566. x0[0] = 1.
  567. assert_equal(vf.fun(x0), ex.fun(x0))
  568. assert x0 is not vf.x
  569. # now test with a HessianUpdate strategy specified
  570. hess = BFGS()
  571. x0 = np.array([1.0, 0.0])
  572. vf = VectorFunction(ex.fun, x0, '3-point', hess, None, None,
  573. (-np.inf, np.inf), None)
  574. with pytest.warns(UserWarning):
  575. # filter UserWarning because ExVectorialFunction is linear and
  576. # a quasi-Newton approximation is used for the Hessian.
  577. assert x0 is not vf.x
  578. assert_equal(vf.fun(x0), ex.fun(x0))
  579. assert x0 is not vf.x
  580. x0[0] = 2.
  581. assert_equal(vf.fun(x0), ex.fun(x0))
  582. assert x0 is not vf.x
  583. x0[0] = 1.
  584. assert_equal(vf.fun(x0), ex.fun(x0))
  585. assert x0 is not vf.x
  586. def test_LinearVectorFunction():
  587. A_dense = np.array([
  588. [-1, 2, 0],
  589. [0, 4, 2]
  590. ])
  591. x0 = np.zeros(3)
  592. A_sparse = csr_matrix(A_dense)
  593. x = np.array([1, -1, 0])
  594. v = np.array([-1, 1])
  595. Ax = np.array([-3, -4])
  596. f1 = LinearVectorFunction(A_dense, x0, None)
  597. assert_(not f1.sparse_jacobian)
  598. f2 = LinearVectorFunction(A_dense, x0, True)
  599. assert_(f2.sparse_jacobian)
  600. f3 = LinearVectorFunction(A_dense, x0, False)
  601. assert_(not f3.sparse_jacobian)
  602. f4 = LinearVectorFunction(A_sparse, x0, None)
  603. assert_(f4.sparse_jacobian)
  604. f5 = LinearVectorFunction(A_sparse, x0, True)
  605. assert_(f5.sparse_jacobian)
  606. f6 = LinearVectorFunction(A_sparse, x0, False)
  607. assert_(not f6.sparse_jacobian)
  608. assert_array_equal(f1.fun(x), Ax)
  609. assert_array_equal(f2.fun(x), Ax)
  610. assert_array_equal(f1.jac(x), A_dense)
  611. assert_array_equal(f2.jac(x).toarray(), A_sparse.toarray())
  612. assert_array_equal(f1.hess(x, v).toarray(), np.zeros((3, 3)))
  613. def test_LinearVectorFunction_memoization():
  614. A = np.array([[-1, 2, 0], [0, 4, 2]])
  615. x0 = np.array([1, 2, -1])
  616. fun = LinearVectorFunction(A, x0, False)
  617. assert_array_equal(x0, fun.x)
  618. assert_array_equal(A.dot(x0), fun.f)
  619. x1 = np.array([-1, 3, 10])
  620. assert_array_equal(A, fun.jac(x1))
  621. assert_array_equal(x1, fun.x)
  622. assert_array_equal(A.dot(x0), fun.f)
  623. assert_array_equal(A.dot(x1), fun.fun(x1))
  624. assert_array_equal(A.dot(x1), fun.f)
  625. def test_IdentityVectorFunction():
  626. x0 = np.zeros(3)
  627. f1 = IdentityVectorFunction(x0, None)
  628. f2 = IdentityVectorFunction(x0, False)
  629. f3 = IdentityVectorFunction(x0, True)
  630. assert_(f1.sparse_jacobian)
  631. assert_(not f2.sparse_jacobian)
  632. assert_(f3.sparse_jacobian)
  633. x = np.array([-1, 2, 1])
  634. v = np.array([-2, 3, 0])
  635. assert_array_equal(f1.fun(x), x)
  636. assert_array_equal(f2.fun(x), x)
  637. assert_array_equal(f1.jac(x).toarray(), np.eye(3))
  638. assert_array_equal(f2.jac(x), np.eye(3))
  639. assert_array_equal(f1.hess(x, v).toarray(), np.zeros((3, 3)))