test_ivp.py 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818
  1. from __future__ import division, print_function, absolute_import
  2. from itertools import product
  3. from numpy.testing import (assert_, assert_allclose,
  4. assert_equal, assert_no_warnings)
  5. import pytest
  6. from pytest import raises as assert_raises
  7. from scipy._lib._numpy_compat import suppress_warnings
  8. import numpy as np
  9. from scipy.optimize._numdiff import group_columns
  10. from scipy.integrate import solve_ivp, RK23, RK45, Radau, BDF, LSODA
  11. from scipy.integrate import OdeSolution
  12. from scipy.integrate._ivp.common import num_jac
  13. from scipy.integrate._ivp.base import ConstantDenseOutput
  14. from scipy.sparse import coo_matrix, csc_matrix
  15. def fun_linear(t, y):
  16. return np.array([-y[0] - 5 * y[1], y[0] + y[1]])
  17. def jac_linear():
  18. return np.array([[-1, -5], [1, 1]])
  19. def sol_linear(t):
  20. return np.vstack((-5 * np.sin(2 * t),
  21. 2 * np.cos(2 * t) + np.sin(2 * t)))
  22. def fun_rational(t, y):
  23. return np.array([y[1] / t,
  24. y[1] * (y[0] + 2 * y[1] - 1) / (t * (y[0] - 1))])
  25. def fun_rational_vectorized(t, y):
  26. return np.vstack((y[1] / t,
  27. y[1] * (y[0] + 2 * y[1] - 1) / (t * (y[0] - 1))))
  28. def jac_rational(t, y):
  29. return np.array([
  30. [0, 1 / t],
  31. [-2 * y[1] ** 2 / (t * (y[0] - 1) ** 2),
  32. (y[0] + 4 * y[1] - 1) / (t * (y[0] - 1))]
  33. ])
  34. def jac_rational_sparse(t, y):
  35. return csc_matrix([
  36. [0, 1 / t],
  37. [-2 * y[1] ** 2 / (t * (y[0] - 1) ** 2),
  38. (y[0] + 4 * y[1] - 1) / (t * (y[0] - 1))]
  39. ])
  40. def sol_rational(t):
  41. return np.asarray((t / (t + 10), 10 * t / (t + 10) ** 2))
  42. def fun_medazko(t, y):
  43. n = y.shape[0] // 2
  44. k = 100
  45. c = 4
  46. phi = 2 if t <= 5 else 0
  47. y = np.hstack((phi, 0, y, y[-2]))
  48. d = 1 / n
  49. j = np.arange(n) + 1
  50. alpha = 2 * (j * d - 1) ** 3 / c ** 2
  51. beta = (j * d - 1) ** 4 / c ** 2
  52. j_2_p1 = 2 * j + 2
  53. j_2_m3 = 2 * j - 2
  54. j_2_m1 = 2 * j
  55. j_2 = 2 * j + 1
  56. f = np.empty(2 * n)
  57. f[::2] = (alpha * (y[j_2_p1] - y[j_2_m3]) / (2 * d) +
  58. beta * (y[j_2_m3] - 2 * y[j_2_m1] + y[j_2_p1]) / d ** 2 -
  59. k * y[j_2_m1] * y[j_2])
  60. f[1::2] = -k * y[j_2] * y[j_2_m1]
  61. return f
  62. def medazko_sparsity(n):
  63. cols = []
  64. rows = []
  65. i = np.arange(n) * 2
  66. cols.append(i[1:])
  67. rows.append(i[1:] - 2)
  68. cols.append(i)
  69. rows.append(i)
  70. cols.append(i)
  71. rows.append(i + 1)
  72. cols.append(i[:-1])
  73. rows.append(i[:-1] + 2)
  74. i = np.arange(n) * 2 + 1
  75. cols.append(i)
  76. rows.append(i)
  77. cols.append(i)
  78. rows.append(i - 1)
  79. cols = np.hstack(cols)
  80. rows = np.hstack(rows)
  81. return coo_matrix((np.ones_like(cols), (cols, rows)))
  82. def fun_complex(t, y):
  83. return -y
  84. def jac_complex(t, y):
  85. return -np.eye(y.shape[0])
  86. def jac_complex_sparse(t, y):
  87. return csc_matrix(jac_complex(t, y))
  88. def sol_complex(t):
  89. y = (0.5 + 1j) * np.exp(-t)
  90. return y.reshape((1, -1))
  91. def compute_error(y, y_true, rtol, atol):
  92. e = (y - y_true) / (atol + rtol * np.abs(y_true))
  93. return np.sqrt(np.sum(np.real(e * e.conj()), axis=0) / e.shape[0])
  94. def test_integration():
  95. rtol = 1e-3
  96. atol = 1e-6
  97. y0 = [1/3, 2/9]
  98. for vectorized, method, t_span, jac in product(
  99. [False, True],
  100. ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA'],
  101. [[5, 9], [5, 1]],
  102. [None, jac_rational, jac_rational_sparse]):
  103. if vectorized:
  104. fun = fun_rational_vectorized
  105. else:
  106. fun = fun_rational
  107. with suppress_warnings() as sup:
  108. sup.filter(UserWarning,
  109. "The following arguments have no effect for a chosen solver: `jac`")
  110. res = solve_ivp(fun, t_span, y0, rtol=rtol,
  111. atol=atol, method=method, dense_output=True,
  112. jac=jac, vectorized=vectorized)
  113. assert_equal(res.t[0], t_span[0])
  114. assert_(res.t_events is None)
  115. assert_(res.success)
  116. assert_equal(res.status, 0)
  117. assert_(res.nfev < 40)
  118. if method in ['RK23', 'RK45', 'LSODA']:
  119. assert_equal(res.njev, 0)
  120. assert_equal(res.nlu, 0)
  121. else:
  122. assert_(0 < res.njev < 3)
  123. assert_(0 < res.nlu < 10)
  124. y_true = sol_rational(res.t)
  125. e = compute_error(res.y, y_true, rtol, atol)
  126. assert_(np.all(e < 5))
  127. tc = np.linspace(*t_span)
  128. yc_true = sol_rational(tc)
  129. yc = res.sol(tc)
  130. e = compute_error(yc, yc_true, rtol, atol)
  131. assert_(np.all(e < 5))
  132. tc = (t_span[0] + t_span[-1]) / 2
  133. yc_true = sol_rational(tc)
  134. yc = res.sol(tc)
  135. e = compute_error(yc, yc_true, rtol, atol)
  136. assert_(np.all(e < 5))
  137. # LSODA for some reasons doesn't pass the polynomial through the
  138. # previous points exactly after the order change. It might be some
  139. # bug in LSOSA implementation or maybe we missing something.
  140. if method != 'LSODA':
  141. assert_allclose(res.sol(res.t), res.y, rtol=1e-15, atol=1e-15)
  142. def test_integration_complex():
  143. rtol = 1e-3
  144. atol = 1e-6
  145. y0 = [0.5 + 1j]
  146. t_span = [0, 1]
  147. tc = np.linspace(t_span[0], t_span[1])
  148. for method, jac in product(['RK23', 'RK45', 'BDF'],
  149. [None, jac_complex, jac_complex_sparse]):
  150. with suppress_warnings() as sup:
  151. sup.filter(UserWarning,
  152. "The following arguments have no effect for a chosen solver: `jac`")
  153. res = solve_ivp(fun_complex, t_span, y0, method=method,
  154. dense_output=True, rtol=rtol, atol=atol, jac=jac)
  155. assert_equal(res.t[0], t_span[0])
  156. assert_(res.t_events is None)
  157. assert_(res.success)
  158. assert_equal(res.status, 0)
  159. assert_(res.nfev < 25)
  160. if method == 'BDF':
  161. assert_equal(res.njev, 1)
  162. assert_(res.nlu < 6)
  163. else:
  164. assert_equal(res.njev, 0)
  165. assert_equal(res.nlu, 0)
  166. y_true = sol_complex(res.t)
  167. e = compute_error(res.y, y_true, rtol, atol)
  168. assert_(np.all(e < 5))
  169. yc_true = sol_complex(tc)
  170. yc = res.sol(tc)
  171. e = compute_error(yc, yc_true, rtol, atol)
  172. assert_(np.all(e < 5))
  173. def test_integration_sparse_difference():
  174. n = 200
  175. t_span = [0, 20]
  176. y0 = np.zeros(2 * n)
  177. y0[1::2] = 1
  178. sparsity = medazko_sparsity(n)
  179. for method in ['BDF', 'Radau']:
  180. res = solve_ivp(fun_medazko, t_span, y0, method=method,
  181. jac_sparsity=sparsity)
  182. assert_equal(res.t[0], t_span[0])
  183. assert_(res.t_events is None)
  184. assert_(res.success)
  185. assert_equal(res.status, 0)
  186. assert_allclose(res.y[78, -1], 0.233994e-3, rtol=1e-2)
  187. assert_allclose(res.y[79, -1], 0, atol=1e-3)
  188. assert_allclose(res.y[148, -1], 0.359561e-3, rtol=1e-2)
  189. assert_allclose(res.y[149, -1], 0, atol=1e-3)
  190. assert_allclose(res.y[198, -1], 0.117374129e-3, rtol=1e-2)
  191. assert_allclose(res.y[199, -1], 0.6190807e-5, atol=1e-3)
  192. assert_allclose(res.y[238, -1], 0, atol=1e-3)
  193. assert_allclose(res.y[239, -1], 0.9999997, rtol=1e-2)
  194. def test_integration_const_jac():
  195. rtol = 1e-3
  196. atol = 1e-6
  197. y0 = [0, 2]
  198. t_span = [0, 2]
  199. J = jac_linear()
  200. J_sparse = csc_matrix(J)
  201. for method, jac in product(['Radau', 'BDF'], [J, J_sparse]):
  202. res = solve_ivp(fun_linear, t_span, y0, rtol=rtol, atol=atol,
  203. method=method, dense_output=True, jac=jac)
  204. assert_equal(res.t[0], t_span[0])
  205. assert_(res.t_events is None)
  206. assert_(res.success)
  207. assert_equal(res.status, 0)
  208. assert_(res.nfev < 100)
  209. assert_equal(res.njev, 0)
  210. assert_(0 < res.nlu < 15)
  211. y_true = sol_linear(res.t)
  212. e = compute_error(res.y, y_true, rtol, atol)
  213. assert_(np.all(e < 10))
  214. tc = np.linspace(*t_span)
  215. yc_true = sol_linear(tc)
  216. yc = res.sol(tc)
  217. e = compute_error(yc, yc_true, rtol, atol)
  218. assert_(np.all(e < 15))
  219. assert_allclose(res.sol(res.t), res.y, rtol=1e-14, atol=1e-14)
  220. @pytest.mark.slow
  221. @pytest.mark.parametrize('method', ['Radau', 'BDF', 'LSODA'])
  222. def test_integration_stiff(method):
  223. rtol = 1e-6
  224. atol = 1e-6
  225. y0 = [1e4, 0, 0]
  226. tspan = [0, 1e8]
  227. def fun_robertson(t, state):
  228. x, y, z = state
  229. return [
  230. -0.04 * x + 1e4 * y * z,
  231. 0.04 * x - 1e4 * y * z - 3e7 * y * y,
  232. 3e7 * y * y,
  233. ]
  234. res = solve_ivp(fun_robertson, tspan, y0, rtol=rtol,
  235. atol=atol, method=method)
  236. # If the stiff mode is not activated correctly, these numbers will be much bigger
  237. assert res.nfev < 5000
  238. assert res.njev < 200
  239. def test_events():
  240. def event_rational_1(t, y):
  241. return y[0] - y[1] ** 0.7
  242. def event_rational_2(t, y):
  243. return y[1] ** 0.6 - y[0]
  244. def event_rational_3(t, y):
  245. return t - 7.4
  246. event_rational_3.terminal = True
  247. for method in ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA']:
  248. res = solve_ivp(fun_rational, [5, 8], [1/3, 2/9], method=method,
  249. events=(event_rational_1, event_rational_2))
  250. assert_equal(res.status, 0)
  251. assert_equal(res.t_events[0].size, 1)
  252. assert_equal(res.t_events[1].size, 1)
  253. assert_(5.3 < res.t_events[0][0] < 5.7)
  254. assert_(7.3 < res.t_events[1][0] < 7.7)
  255. event_rational_1.direction = 1
  256. event_rational_2.direction = 1
  257. res = solve_ivp(fun_rational, [5, 8], [1 / 3, 2 / 9], method=method,
  258. events=(event_rational_1, event_rational_2))
  259. assert_equal(res.status, 0)
  260. assert_equal(res.t_events[0].size, 1)
  261. assert_equal(res.t_events[1].size, 0)
  262. assert_(5.3 < res.t_events[0][0] < 5.7)
  263. event_rational_1.direction = -1
  264. event_rational_2.direction = -1
  265. res = solve_ivp(fun_rational, [5, 8], [1 / 3, 2 / 9], method=method,
  266. events=(event_rational_1, event_rational_2))
  267. assert_equal(res.status, 0)
  268. assert_equal(res.t_events[0].size, 0)
  269. assert_equal(res.t_events[1].size, 1)
  270. assert_(7.3 < res.t_events[1][0] < 7.7)
  271. event_rational_1.direction = 0
  272. event_rational_2.direction = 0
  273. res = solve_ivp(fun_rational, [5, 8], [1 / 3, 2 / 9], method=method,
  274. events=(event_rational_1, event_rational_2,
  275. event_rational_3), dense_output=True)
  276. assert_equal(res.status, 1)
  277. assert_equal(res.t_events[0].size, 1)
  278. assert_equal(res.t_events[1].size, 0)
  279. assert_equal(res.t_events[2].size, 1)
  280. assert_(5.3 < res.t_events[0][0] < 5.7)
  281. assert_(7.3 < res.t_events[2][0] < 7.5)
  282. res = solve_ivp(fun_rational, [5, 8], [1 / 3, 2 / 9], method=method,
  283. events=event_rational_1, dense_output=True)
  284. assert_equal(res.status, 0)
  285. assert_equal(res.t_events[0].size, 1)
  286. assert_(5.3 < res.t_events[0][0] < 5.7)
  287. # Also test that termination by event doesn't break interpolants.
  288. tc = np.linspace(res.t[0], res.t[-1])
  289. yc_true = sol_rational(tc)
  290. yc = res.sol(tc)
  291. e = compute_error(yc, yc_true, 1e-3, 1e-6)
  292. assert_(np.all(e < 5))
  293. # Test in backward direction.
  294. event_rational_1.direction = 0
  295. event_rational_2.direction = 0
  296. for method in ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA']:
  297. res = solve_ivp(fun_rational, [8, 5], [4/9, 20/81], method=method,
  298. events=(event_rational_1, event_rational_2))
  299. assert_equal(res.status, 0)
  300. assert_equal(res.t_events[0].size, 1)
  301. assert_equal(res.t_events[1].size, 1)
  302. assert_(5.3 < res.t_events[0][0] < 5.7)
  303. assert_(7.3 < res.t_events[1][0] < 7.7)
  304. event_rational_1.direction = -1
  305. event_rational_2.direction = -1
  306. res = solve_ivp(fun_rational, [8, 5], [4/9, 20/81], method=method,
  307. events=(event_rational_1, event_rational_2))
  308. assert_equal(res.status, 0)
  309. assert_equal(res.t_events[0].size, 1)
  310. assert_equal(res.t_events[1].size, 0)
  311. assert_(5.3 < res.t_events[0][0] < 5.7)
  312. event_rational_1.direction = 1
  313. event_rational_2.direction = 1
  314. res = solve_ivp(fun_rational, [8, 5], [4/9, 20/81], method=method,
  315. events=(event_rational_1, event_rational_2))
  316. assert_equal(res.status, 0)
  317. assert_equal(res.t_events[0].size, 0)
  318. assert_equal(res.t_events[1].size, 1)
  319. assert_(7.3 < res.t_events[1][0] < 7.7)
  320. event_rational_1.direction = 0
  321. event_rational_2.direction = 0
  322. res = solve_ivp(fun_rational, [8, 5], [4/9, 20/81], method=method,
  323. events=(event_rational_1, event_rational_2,
  324. event_rational_3), dense_output=True)
  325. assert_equal(res.status, 1)
  326. assert_equal(res.t_events[0].size, 0)
  327. assert_equal(res.t_events[1].size, 1)
  328. assert_equal(res.t_events[2].size, 1)
  329. assert_(7.3 < res.t_events[1][0] < 7.7)
  330. assert_(7.3 < res.t_events[2][0] < 7.5)
  331. # Also test that termination by event doesn't break interpolants.
  332. tc = np.linspace(res.t[-1], res.t[0])
  333. yc_true = sol_rational(tc)
  334. yc = res.sol(tc)
  335. e = compute_error(yc, yc_true, 1e-3, 1e-6)
  336. assert_(np.all(e < 5))
  337. def test_max_step():
  338. rtol = 1e-3
  339. atol = 1e-6
  340. y0 = [1/3, 2/9]
  341. for method in [RK23, RK45, Radau, BDF, LSODA]:
  342. for t_span in ([5, 9], [5, 1]):
  343. res = solve_ivp(fun_rational, t_span, y0, rtol=rtol,
  344. max_step=0.5, atol=atol, method=method,
  345. dense_output=True)
  346. assert_equal(res.t[0], t_span[0])
  347. assert_equal(res.t[-1], t_span[-1])
  348. assert_(np.all(np.abs(np.diff(res.t)) <= 0.5))
  349. assert_(res.t_events is None)
  350. assert_(res.success)
  351. assert_equal(res.status, 0)
  352. y_true = sol_rational(res.t)
  353. e = compute_error(res.y, y_true, rtol, atol)
  354. assert_(np.all(e < 5))
  355. tc = np.linspace(*t_span)
  356. yc_true = sol_rational(tc)
  357. yc = res.sol(tc)
  358. e = compute_error(yc, yc_true, rtol, atol)
  359. assert_(np.all(e < 5))
  360. # See comment in test_integration.
  361. if method is not LSODA:
  362. assert_allclose(res.sol(res.t), res.y, rtol=1e-15, atol=1e-15)
  363. assert_raises(ValueError, method, fun_rational, t_span[0], y0,
  364. t_span[1], max_step=-1)
  365. if method is not LSODA:
  366. solver = method(fun_rational, t_span[0], y0, t_span[1],
  367. rtol=rtol, atol=atol, max_step=1e-20)
  368. message = solver.step()
  369. assert_equal(solver.status, 'failed')
  370. assert_("step size is less" in message)
  371. assert_raises(RuntimeError, solver.step)
  372. def test_first_step():
  373. rtol = 1e-3
  374. atol = 1e-6
  375. y0 = [1/3, 2/9]
  376. first_step = 0.1
  377. for method in [RK23, RK45, Radau, BDF, LSODA]:
  378. for t_span in ([5, 9], [5, 1]):
  379. res = solve_ivp(fun_rational, t_span, y0, rtol=rtol,
  380. max_step=0.5, atol=atol, method=method,
  381. dense_output=True, first_step=first_step)
  382. assert_equal(res.t[0], t_span[0])
  383. assert_equal(res.t[-1], t_span[-1])
  384. assert_allclose(first_step, np.abs(res.t[1] - 5))
  385. assert_(res.t_events is None)
  386. assert_(res.success)
  387. assert_equal(res.status, 0)
  388. y_true = sol_rational(res.t)
  389. e = compute_error(res.y, y_true, rtol, atol)
  390. assert_(np.all(e < 5))
  391. tc = np.linspace(*t_span)
  392. yc_true = sol_rational(tc)
  393. yc = res.sol(tc)
  394. e = compute_error(yc, yc_true, rtol, atol)
  395. assert_(np.all(e < 5))
  396. # See comment in test_integration.
  397. if method is not LSODA:
  398. assert_allclose(res.sol(res.t), res.y, rtol=1e-15, atol=1e-15)
  399. assert_raises(ValueError, method, fun_rational, t_span[0], y0,
  400. t_span[1], first_step=-1)
  401. assert_raises(ValueError, method, fun_rational, t_span[0], y0,
  402. t_span[1], first_step=5)
  403. def test_t_eval():
  404. rtol = 1e-3
  405. atol = 1e-6
  406. y0 = [1/3, 2/9]
  407. for t_span in ([5, 9], [5, 1]):
  408. t_eval = np.linspace(t_span[0], t_span[1], 10)
  409. res = solve_ivp(fun_rational, t_span, y0, rtol=rtol, atol=atol,
  410. t_eval=t_eval)
  411. assert_equal(res.t, t_eval)
  412. assert_(res.t_events is None)
  413. assert_(res.success)
  414. assert_equal(res.status, 0)
  415. y_true = sol_rational(res.t)
  416. e = compute_error(res.y, y_true, rtol, atol)
  417. assert_(np.all(e < 5))
  418. t_eval = [5, 5.01, 7, 8, 8.01, 9]
  419. res = solve_ivp(fun_rational, [5, 9], y0, rtol=rtol, atol=atol,
  420. t_eval=t_eval)
  421. assert_equal(res.t, t_eval)
  422. assert_(res.t_events is None)
  423. assert_(res.success)
  424. assert_equal(res.status, 0)
  425. y_true = sol_rational(res.t)
  426. e = compute_error(res.y, y_true, rtol, atol)
  427. assert_(np.all(e < 5))
  428. t_eval = [5, 4.99, 3, 1.5, 1.1, 1.01, 1]
  429. res = solve_ivp(fun_rational, [5, 1], y0, rtol=rtol, atol=atol,
  430. t_eval=t_eval)
  431. assert_equal(res.t, t_eval)
  432. assert_(res.t_events is None)
  433. assert_(res.success)
  434. assert_equal(res.status, 0)
  435. t_eval = [5.01, 7, 8, 8.01]
  436. res = solve_ivp(fun_rational, [5, 9], y0, rtol=rtol, atol=atol,
  437. t_eval=t_eval)
  438. assert_equal(res.t, t_eval)
  439. assert_(res.t_events is None)
  440. assert_(res.success)
  441. assert_equal(res.status, 0)
  442. y_true = sol_rational(res.t)
  443. e = compute_error(res.y, y_true, rtol, atol)
  444. assert_(np.all(e < 5))
  445. t_eval = [4.99, 3, 1.5, 1.1, 1.01]
  446. res = solve_ivp(fun_rational, [5, 1], y0, rtol=rtol, atol=atol,
  447. t_eval=t_eval)
  448. assert_equal(res.t, t_eval)
  449. assert_(res.t_events is None)
  450. assert_(res.success)
  451. assert_equal(res.status, 0)
  452. t_eval = [4, 6]
  453. assert_raises(ValueError, solve_ivp, fun_rational, [5, 9], y0,
  454. rtol=rtol, atol=atol, t_eval=t_eval)
  455. def test_t_eval_dense_output():
  456. rtol = 1e-3
  457. atol = 1e-6
  458. y0 = [1/3, 2/9]
  459. t_span = [5, 9]
  460. t_eval = np.linspace(t_span[0], t_span[1], 10)
  461. res = solve_ivp(fun_rational, t_span, y0, rtol=rtol, atol=atol,
  462. t_eval=t_eval)
  463. res_d = solve_ivp(fun_rational, t_span, y0, rtol=rtol, atol=atol,
  464. t_eval=t_eval, dense_output=True)
  465. assert_equal(res.t, t_eval)
  466. assert_(res.t_events is None)
  467. assert_(res.success)
  468. assert_equal(res.status, 0)
  469. assert_equal(res.t, res_d.t)
  470. assert_equal(res.y, res_d.y)
  471. assert_(res_d.t_events is None)
  472. assert_(res_d.success)
  473. assert_equal(res_d.status, 0)
  474. # if t and y are equal only test values for one case
  475. y_true = sol_rational(res.t)
  476. e = compute_error(res.y, y_true, rtol, atol)
  477. assert_(np.all(e < 5))
  478. def test_no_integration():
  479. for method in ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA']:
  480. sol = solve_ivp(lambda t, y: -y, [4, 4], [2, 3],
  481. method=method, dense_output=True)
  482. assert_equal(sol.sol(4), [2, 3])
  483. assert_equal(sol.sol([4, 5, 6]), [[2, 2, 2], [3, 3, 3]])
  484. def test_no_integration_class():
  485. for method in [RK23, RK45, Radau, BDF, LSODA]:
  486. solver = method(lambda t, y: -y, 0.0, [10.0, 0.0], 0.0)
  487. solver.step()
  488. assert_equal(solver.status, 'finished')
  489. sol = solver.dense_output()
  490. assert_equal(sol(0.0), [10.0, 0.0])
  491. assert_equal(sol([0, 1, 2]), [[10, 10, 10], [0, 0, 0]])
  492. solver = method(lambda t, y: -y, 0.0, [], np.inf)
  493. solver.step()
  494. assert_equal(solver.status, 'finished')
  495. sol = solver.dense_output()
  496. assert_equal(sol(100.0), [])
  497. assert_equal(sol([0, 1, 2]), np.empty((0, 3)))
  498. def test_empty():
  499. def fun(t, y):
  500. return np.zeros((0,))
  501. y0 = np.zeros((0,))
  502. for method in ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA']:
  503. sol = assert_no_warnings(solve_ivp, fun, [0, 10], y0,
  504. method=method, dense_output=True)
  505. assert_equal(sol.sol(10), np.zeros((0,)))
  506. assert_equal(sol.sol([1, 2, 3]), np.zeros((0, 3)))
  507. for method in ['RK23', 'RK45', 'Radau', 'BDF', 'LSODA']:
  508. sol = assert_no_warnings(solve_ivp, fun, [0, np.inf], y0,
  509. method=method, dense_output=True)
  510. assert_equal(sol.sol(10), np.zeros((0,)))
  511. assert_equal(sol.sol([1, 2, 3]), np.zeros((0, 3)))
  512. def test_ConstantDenseOutput():
  513. sol = ConstantDenseOutput(0, 1, np.array([1, 2]))
  514. assert_allclose(sol(1.5), [1, 2])
  515. assert_allclose(sol([1, 1.5, 2]), [[1, 1, 1], [2, 2, 2]])
  516. sol = ConstantDenseOutput(0, 1, np.array([]))
  517. assert_allclose(sol(1.5), np.empty(0))
  518. assert_allclose(sol([1, 1.5, 2]), np.empty((0, 3)))
  519. def test_classes():
  520. y0 = [1 / 3, 2 / 9]
  521. for cls in [RK23, RK45, Radau, BDF, LSODA]:
  522. solver = cls(fun_rational, 5, y0, np.inf)
  523. assert_equal(solver.n, 2)
  524. assert_equal(solver.status, 'running')
  525. assert_equal(solver.t_bound, np.inf)
  526. assert_equal(solver.direction, 1)
  527. assert_equal(solver.t, 5)
  528. assert_equal(solver.y, y0)
  529. assert_(solver.step_size is None)
  530. if cls is not LSODA:
  531. assert_(solver.nfev > 0)
  532. assert_(solver.njev >= 0)
  533. assert_equal(solver.nlu, 0)
  534. else:
  535. assert_equal(solver.nfev, 0)
  536. assert_equal(solver.njev, 0)
  537. assert_equal(solver.nlu, 0)
  538. assert_raises(RuntimeError, solver.dense_output)
  539. message = solver.step()
  540. assert_equal(solver.status, 'running')
  541. assert_equal(message, None)
  542. assert_equal(solver.n, 2)
  543. assert_equal(solver.t_bound, np.inf)
  544. assert_equal(solver.direction, 1)
  545. assert_(solver.t > 5)
  546. assert_(not np.all(np.equal(solver.y, y0)))
  547. assert_(solver.step_size > 0)
  548. assert_(solver.nfev > 0)
  549. assert_(solver.njev >= 0)
  550. assert_(solver.nlu >= 0)
  551. sol = solver.dense_output()
  552. assert_allclose(sol(5), y0, rtol=1e-15, atol=0)
  553. def test_OdeSolution():
  554. ts = np.array([0, 2, 5], dtype=float)
  555. s1 = ConstantDenseOutput(ts[0], ts[1], np.array([-1]))
  556. s2 = ConstantDenseOutput(ts[1], ts[2], np.array([1]))
  557. sol = OdeSolution(ts, [s1, s2])
  558. assert_equal(sol(-1), [-1])
  559. assert_equal(sol(1), [-1])
  560. assert_equal(sol(2), [-1])
  561. assert_equal(sol(3), [1])
  562. assert_equal(sol(5), [1])
  563. assert_equal(sol(6), [1])
  564. assert_equal(sol([0, 6, -2, 1.5, 4.5, 2.5, 5, 5.5, 2]),
  565. np.array([[-1, 1, -1, -1, 1, 1, 1, 1, -1]]))
  566. ts = np.array([10, 4, -3])
  567. s1 = ConstantDenseOutput(ts[0], ts[1], np.array([-1]))
  568. s2 = ConstantDenseOutput(ts[1], ts[2], np.array([1]))
  569. sol = OdeSolution(ts, [s1, s2])
  570. assert_equal(sol(11), [-1])
  571. assert_equal(sol(10), [-1])
  572. assert_equal(sol(5), [-1])
  573. assert_equal(sol(4), [-1])
  574. assert_equal(sol(0), [1])
  575. assert_equal(sol(-3), [1])
  576. assert_equal(sol(-4), [1])
  577. assert_equal(sol([12, -5, 10, -3, 6, 1, 4]),
  578. np.array([[-1, 1, -1, 1, -1, 1, -1]]))
  579. ts = np.array([1, 1])
  580. s = ConstantDenseOutput(1, 1, np.array([10]))
  581. sol = OdeSolution(ts, [s])
  582. assert_equal(sol(0), [10])
  583. assert_equal(sol(1), [10])
  584. assert_equal(sol(2), [10])
  585. assert_equal(sol([2, 1, 0]), np.array([[10, 10, 10]]))
  586. def test_num_jac():
  587. def fun(t, y):
  588. return np.vstack([
  589. -0.04 * y[0] + 1e4 * y[1] * y[2],
  590. 0.04 * y[0] - 1e4 * y[1] * y[2] - 3e7 * y[1] ** 2,
  591. 3e7 * y[1] ** 2
  592. ])
  593. def jac(t, y):
  594. return np.array([
  595. [-0.04, 1e4 * y[2], 1e4 * y[1]],
  596. [0.04, -1e4 * y[2] - 6e7 * y[1], -1e4 * y[1]],
  597. [0, 6e7 * y[1], 0]
  598. ])
  599. t = 1
  600. y = np.array([1, 0, 0])
  601. J_true = jac(t, y)
  602. threshold = 1e-5
  603. f = fun(t, y).ravel()
  604. J_num, factor = num_jac(fun, t, y, f, threshold, None)
  605. assert_allclose(J_num, J_true, rtol=1e-5, atol=1e-5)
  606. J_num, factor = num_jac(fun, t, y, f, threshold, factor)
  607. assert_allclose(J_num, J_true, rtol=1e-5, atol=1e-5)
  608. def test_num_jac_sparse():
  609. def fun(t, y):
  610. e = y[1:]**3 - y[:-1]**2
  611. z = np.zeros(y.shape[1])
  612. return np.vstack((z, 3 * e)) + np.vstack((2 * e, z))
  613. def structure(n):
  614. A = np.zeros((n, n), dtype=int)
  615. A[0, 0] = 1
  616. A[0, 1] = 1
  617. for i in range(1, n - 1):
  618. A[i, i - 1: i + 2] = 1
  619. A[-1, -1] = 1
  620. A[-1, -2] = 1
  621. return A
  622. np.random.seed(0)
  623. n = 20
  624. y = np.random.randn(n)
  625. A = structure(n)
  626. groups = group_columns(A)
  627. f = fun(0, y[:, None]).ravel()
  628. # Compare dense and sparse results, assuming that dense implementation
  629. # is correct (as it is straightforward).
  630. J_num_sparse, factor_sparse = num_jac(fun, 0, y.ravel(), f, 1e-8, None,
  631. sparsity=(A, groups))
  632. J_num_dense, factor_dense = num_jac(fun, 0, y.ravel(), f, 1e-8, None)
  633. assert_allclose(J_num_dense, J_num_sparse.toarray(),
  634. rtol=1e-12, atol=1e-14)
  635. assert_allclose(factor_dense, factor_sparse, rtol=1e-12, atol=1e-14)
  636. # Take small factors to trigger their recomputing inside.
  637. factor = np.random.uniform(0, 1e-12, size=n)
  638. J_num_sparse, factor_sparse = num_jac(fun, 0, y.ravel(), f, 1e-8, factor,
  639. sparsity=(A, groups))
  640. J_num_dense, factor_dense = num_jac(fun, 0, y.ravel(), f, 1e-8, factor)
  641. assert_allclose(J_num_dense, J_num_sparse.toarray(),
  642. rtol=1e-12, atol=1e-14)
  643. assert_allclose(factor_dense, factor_sparse, rtol=1e-12, atol=1e-14)