test_integrate.py 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830
  1. # Authors: Nils Wagner, Ed Schofield, Pauli Virtanen, John Travers
  2. """
  3. Tests for numerical integration.
  4. """
  5. import numpy as np
  6. from numpy import (arange, zeros, array, dot, sqrt, cos, sin, eye, pi, exp,
  7. allclose)
  8. from numpy.testing import (
  9. assert_, assert_array_almost_equal,
  10. assert_allclose, assert_array_equal, assert_equal, assert_warns)
  11. from pytest import raises as assert_raises
  12. from scipy.integrate import odeint, ode, complex_ode
  13. #------------------------------------------------------------------------------
  14. # Test ODE integrators
  15. #------------------------------------------------------------------------------
  16. class TestOdeint:
  17. # Check integrate.odeint
  18. def _do_problem(self, problem):
  19. t = arange(0.0, problem.stop_t, 0.05)
  20. # Basic case
  21. z, infodict = odeint(problem.f, problem.z0, t, full_output=True)
  22. assert_(problem.verify(z, t))
  23. # Use tfirst=True
  24. z, infodict = odeint(lambda t, y: problem.f(y, t), problem.z0, t,
  25. full_output=True, tfirst=True)
  26. assert_(problem.verify(z, t))
  27. if hasattr(problem, 'jac'):
  28. # Use Dfun
  29. z, infodict = odeint(problem.f, problem.z0, t, Dfun=problem.jac,
  30. full_output=True)
  31. assert_(problem.verify(z, t))
  32. # Use Dfun and tfirst=True
  33. z, infodict = odeint(lambda t, y: problem.f(y, t), problem.z0, t,
  34. Dfun=lambda t, y: problem.jac(y, t),
  35. full_output=True, tfirst=True)
  36. assert_(problem.verify(z, t))
  37. def test_odeint(self):
  38. for problem_cls in PROBLEMS:
  39. problem = problem_cls()
  40. if problem.cmplx:
  41. continue
  42. self._do_problem(problem)
  43. class TestODEClass:
  44. ode_class = None # Set in subclass.
  45. def _do_problem(self, problem, integrator, method='adams'):
  46. # ode has callback arguments in different order than odeint
  47. f = lambda t, z: problem.f(z, t)
  48. jac = None
  49. if hasattr(problem, 'jac'):
  50. jac = lambda t, z: problem.jac(z, t)
  51. integrator_params = {}
  52. if problem.lband is not None or problem.uband is not None:
  53. integrator_params['uband'] = problem.uband
  54. integrator_params['lband'] = problem.lband
  55. ig = self.ode_class(f, jac)
  56. ig.set_integrator(integrator,
  57. atol=problem.atol/10,
  58. rtol=problem.rtol/10,
  59. method=method,
  60. **integrator_params)
  61. ig.set_initial_value(problem.z0, t=0.0)
  62. z = ig.integrate(problem.stop_t)
  63. assert_array_equal(z, ig.y)
  64. assert_(ig.successful(), (problem, method))
  65. assert_(ig.get_return_code() > 0, (problem, method))
  66. assert_(problem.verify(array([z]), problem.stop_t), (problem, method))
  67. class TestOde(TestODEClass):
  68. ode_class = ode
  69. def test_vode(self):
  70. # Check the vode solver
  71. for problem_cls in PROBLEMS:
  72. problem = problem_cls()
  73. if problem.cmplx:
  74. continue
  75. if not problem.stiff:
  76. self._do_problem(problem, 'vode', 'adams')
  77. self._do_problem(problem, 'vode', 'bdf')
  78. def test_zvode(self):
  79. # Check the zvode solver
  80. for problem_cls in PROBLEMS:
  81. problem = problem_cls()
  82. if not problem.stiff:
  83. self._do_problem(problem, 'zvode', 'adams')
  84. self._do_problem(problem, 'zvode', 'bdf')
  85. def test_lsoda(self):
  86. # Check the lsoda solver
  87. for problem_cls in PROBLEMS:
  88. problem = problem_cls()
  89. if problem.cmplx:
  90. continue
  91. self._do_problem(problem, 'lsoda')
  92. def test_dopri5(self):
  93. # Check the dopri5 solver
  94. for problem_cls in PROBLEMS:
  95. problem = problem_cls()
  96. if problem.cmplx:
  97. continue
  98. if problem.stiff:
  99. continue
  100. if hasattr(problem, 'jac'):
  101. continue
  102. self._do_problem(problem, 'dopri5')
  103. def test_dop853(self):
  104. # Check the dop853 solver
  105. for problem_cls in PROBLEMS:
  106. problem = problem_cls()
  107. if problem.cmplx:
  108. continue
  109. if problem.stiff:
  110. continue
  111. if hasattr(problem, 'jac'):
  112. continue
  113. self._do_problem(problem, 'dop853')
  114. def test_concurrent_fail(self):
  115. for sol in ('vode', 'zvode', 'lsoda'):
  116. f = lambda t, y: 1.0
  117. r = ode(f).set_integrator(sol)
  118. r.set_initial_value(0, 0)
  119. r2 = ode(f).set_integrator(sol)
  120. r2.set_initial_value(0, 0)
  121. r.integrate(r.t + 0.1)
  122. r2.integrate(r2.t + 0.1)
  123. assert_raises(RuntimeError, r.integrate, r.t + 0.1)
  124. def test_concurrent_ok(self):
  125. f = lambda t, y: 1.0
  126. for k in range(3):
  127. for sol in ('vode', 'zvode', 'lsoda', 'dopri5', 'dop853'):
  128. r = ode(f).set_integrator(sol)
  129. r.set_initial_value(0, 0)
  130. r2 = ode(f).set_integrator(sol)
  131. r2.set_initial_value(0, 0)
  132. r.integrate(r.t + 0.1)
  133. r2.integrate(r2.t + 0.1)
  134. r2.integrate(r2.t + 0.1)
  135. assert_allclose(r.y, 0.1)
  136. assert_allclose(r2.y, 0.2)
  137. for sol in ('dopri5', 'dop853'):
  138. r = ode(f).set_integrator(sol)
  139. r.set_initial_value(0, 0)
  140. r2 = ode(f).set_integrator(sol)
  141. r2.set_initial_value(0, 0)
  142. r.integrate(r.t + 0.1)
  143. r.integrate(r.t + 0.1)
  144. r2.integrate(r2.t + 0.1)
  145. r.integrate(r.t + 0.1)
  146. r2.integrate(r2.t + 0.1)
  147. assert_allclose(r.y, 0.3)
  148. assert_allclose(r2.y, 0.2)
  149. class TestComplexOde(TestODEClass):
  150. ode_class = complex_ode
  151. def test_vode(self):
  152. # Check the vode solver
  153. for problem_cls in PROBLEMS:
  154. problem = problem_cls()
  155. if not problem.stiff:
  156. self._do_problem(problem, 'vode', 'adams')
  157. else:
  158. self._do_problem(problem, 'vode', 'bdf')
  159. def test_lsoda(self):
  160. # Check the lsoda solver
  161. for problem_cls in PROBLEMS:
  162. problem = problem_cls()
  163. self._do_problem(problem, 'lsoda')
  164. def test_dopri5(self):
  165. # Check the dopri5 solver
  166. for problem_cls in PROBLEMS:
  167. problem = problem_cls()
  168. if problem.stiff:
  169. continue
  170. if hasattr(problem, 'jac'):
  171. continue
  172. self._do_problem(problem, 'dopri5')
  173. def test_dop853(self):
  174. # Check the dop853 solver
  175. for problem_cls in PROBLEMS:
  176. problem = problem_cls()
  177. if problem.stiff:
  178. continue
  179. if hasattr(problem, 'jac'):
  180. continue
  181. self._do_problem(problem, 'dop853')
  182. class TestSolout:
  183. # Check integrate.ode correctly handles solout for dopri5 and dop853
  184. def _run_solout_test(self, integrator):
  185. # Check correct usage of solout
  186. ts = []
  187. ys = []
  188. t0 = 0.0
  189. tend = 10.0
  190. y0 = [1.0, 2.0]
  191. def solout(t, y):
  192. ts.append(t)
  193. ys.append(y.copy())
  194. def rhs(t, y):
  195. return [y[0] + y[1], -y[1]**2]
  196. ig = ode(rhs).set_integrator(integrator)
  197. ig.set_solout(solout)
  198. ig.set_initial_value(y0, t0)
  199. ret = ig.integrate(tend)
  200. assert_array_equal(ys[0], y0)
  201. assert_array_equal(ys[-1], ret)
  202. assert_equal(ts[0], t0)
  203. assert_equal(ts[-1], tend)
  204. def test_solout(self):
  205. for integrator in ('dopri5', 'dop853'):
  206. self._run_solout_test(integrator)
  207. def _run_solout_after_initial_test(self, integrator):
  208. # Check if solout works even if it is set after the initial value.
  209. ts = []
  210. ys = []
  211. t0 = 0.0
  212. tend = 10.0
  213. y0 = [1.0, 2.0]
  214. def solout(t, y):
  215. ts.append(t)
  216. ys.append(y.copy())
  217. def rhs(t, y):
  218. return [y[0] + y[1], -y[1]**2]
  219. ig = ode(rhs).set_integrator(integrator)
  220. ig.set_initial_value(y0, t0)
  221. ig.set_solout(solout)
  222. ret = ig.integrate(tend)
  223. assert_array_equal(ys[0], y0)
  224. assert_array_equal(ys[-1], ret)
  225. assert_equal(ts[0], t0)
  226. assert_equal(ts[-1], tend)
  227. def test_solout_after_initial(self):
  228. for integrator in ('dopri5', 'dop853'):
  229. self._run_solout_after_initial_test(integrator)
  230. def _run_solout_break_test(self, integrator):
  231. # Check correct usage of stopping via solout
  232. ts = []
  233. ys = []
  234. t0 = 0.0
  235. tend = 10.0
  236. y0 = [1.0, 2.0]
  237. def solout(t, y):
  238. ts.append(t)
  239. ys.append(y.copy())
  240. if t > tend/2.0:
  241. return -1
  242. def rhs(t, y):
  243. return [y[0] + y[1], -y[1]**2]
  244. ig = ode(rhs).set_integrator(integrator)
  245. ig.set_solout(solout)
  246. ig.set_initial_value(y0, t0)
  247. ret = ig.integrate(tend)
  248. assert_array_equal(ys[0], y0)
  249. assert_array_equal(ys[-1], ret)
  250. assert_equal(ts[0], t0)
  251. assert_(ts[-1] > tend/2.0)
  252. assert_(ts[-1] < tend)
  253. def test_solout_break(self):
  254. for integrator in ('dopri5', 'dop853'):
  255. self._run_solout_break_test(integrator)
  256. class TestComplexSolout:
  257. # Check integrate.ode correctly handles solout for dopri5 and dop853
  258. def _run_solout_test(self, integrator):
  259. # Check correct usage of solout
  260. ts = []
  261. ys = []
  262. t0 = 0.0
  263. tend = 20.0
  264. y0 = [0.0]
  265. def solout(t, y):
  266. ts.append(t)
  267. ys.append(y.copy())
  268. def rhs(t, y):
  269. return [1.0/(t - 10.0 - 1j)]
  270. ig = complex_ode(rhs).set_integrator(integrator)
  271. ig.set_solout(solout)
  272. ig.set_initial_value(y0, t0)
  273. ret = ig.integrate(tend)
  274. assert_array_equal(ys[0], y0)
  275. assert_array_equal(ys[-1], ret)
  276. assert_equal(ts[0], t0)
  277. assert_equal(ts[-1], tend)
  278. def test_solout(self):
  279. for integrator in ('dopri5', 'dop853'):
  280. self._run_solout_test(integrator)
  281. def _run_solout_break_test(self, integrator):
  282. # Check correct usage of stopping via solout
  283. ts = []
  284. ys = []
  285. t0 = 0.0
  286. tend = 20.0
  287. y0 = [0.0]
  288. def solout(t, y):
  289. ts.append(t)
  290. ys.append(y.copy())
  291. if t > tend/2.0:
  292. return -1
  293. def rhs(t, y):
  294. return [1.0/(t - 10.0 - 1j)]
  295. ig = complex_ode(rhs).set_integrator(integrator)
  296. ig.set_solout(solout)
  297. ig.set_initial_value(y0, t0)
  298. ret = ig.integrate(tend)
  299. assert_array_equal(ys[0], y0)
  300. assert_array_equal(ys[-1], ret)
  301. assert_equal(ts[0], t0)
  302. assert_(ts[-1] > tend/2.0)
  303. assert_(ts[-1] < tend)
  304. def test_solout_break(self):
  305. for integrator in ('dopri5', 'dop853'):
  306. self._run_solout_break_test(integrator)
  307. #------------------------------------------------------------------------------
  308. # Test problems
  309. #------------------------------------------------------------------------------
  310. class ODE:
  311. """
  312. ODE problem
  313. """
  314. stiff = False
  315. cmplx = False
  316. stop_t = 1
  317. z0 = []
  318. lband = None
  319. uband = None
  320. atol = 1e-6
  321. rtol = 1e-5
  322. class SimpleOscillator(ODE):
  323. r"""
  324. Free vibration of a simple oscillator::
  325. m \ddot{u} + k u = 0, u(0) = u_0 \dot{u}(0) \dot{u}_0
  326. Solution::
  327. u(t) = u_0*cos(sqrt(k/m)*t)+\dot{u}_0*sin(sqrt(k/m)*t)/sqrt(k/m)
  328. """
  329. stop_t = 1 + 0.09
  330. z0 = array([1.0, 0.1], float)
  331. k = 4.0
  332. m = 1.0
  333. def f(self, z, t):
  334. tmp = zeros((2, 2), float)
  335. tmp[0, 1] = 1.0
  336. tmp[1, 0] = -self.k / self.m
  337. return dot(tmp, z)
  338. def verify(self, zs, t):
  339. omega = sqrt(self.k / self.m)
  340. u = self.z0[0]*cos(omega*t) + self.z0[1]*sin(omega*t)/omega
  341. return allclose(u, zs[:, 0], atol=self.atol, rtol=self.rtol)
  342. class ComplexExp(ODE):
  343. r"""The equation :lm:`\dot u = i u`"""
  344. stop_t = 1.23*pi
  345. z0 = exp([1j, 2j, 3j, 4j, 5j])
  346. cmplx = True
  347. def f(self, z, t):
  348. return 1j*z
  349. def jac(self, z, t):
  350. return 1j*eye(5)
  351. def verify(self, zs, t):
  352. u = self.z0 * exp(1j*t)
  353. return allclose(u, zs, atol=self.atol, rtol=self.rtol)
  354. class Pi(ODE):
  355. r"""Integrate 1/(t + 1j) from t=-10 to t=10"""
  356. stop_t = 20
  357. z0 = [0]
  358. cmplx = True
  359. def f(self, z, t):
  360. return array([1./(t - 10 + 1j)])
  361. def verify(self, zs, t):
  362. u = -2j * np.arctan(10)
  363. return allclose(u, zs[-1, :], atol=self.atol, rtol=self.rtol)
  364. class CoupledDecay(ODE):
  365. r"""
  366. 3 coupled decays suited for banded treatment
  367. (banded mode makes it necessary when N>>3)
  368. """
  369. stiff = True
  370. stop_t = 0.5
  371. z0 = [5.0, 7.0, 13.0]
  372. lband = 1
  373. uband = 0
  374. lmbd = [0.17, 0.23, 0.29] # fictitious decay constants
  375. def f(self, z, t):
  376. lmbd = self.lmbd
  377. return np.array([-lmbd[0]*z[0],
  378. -lmbd[1]*z[1] + lmbd[0]*z[0],
  379. -lmbd[2]*z[2] + lmbd[1]*z[1]])
  380. def jac(self, z, t):
  381. # The full Jacobian is
  382. #
  383. # [-lmbd[0] 0 0 ]
  384. # [ lmbd[0] -lmbd[1] 0 ]
  385. # [ 0 lmbd[1] -lmbd[2]]
  386. #
  387. # The lower and upper bandwidths are lband=1 and uband=0, resp.
  388. # The representation of this array in packed format is
  389. #
  390. # [-lmbd[0] -lmbd[1] -lmbd[2]]
  391. # [ lmbd[0] lmbd[1] 0 ]
  392. lmbd = self.lmbd
  393. j = np.zeros((self.lband + self.uband + 1, 3), order='F')
  394. def set_j(ri, ci, val):
  395. j[self.uband + ri - ci, ci] = val
  396. set_j(0, 0, -lmbd[0])
  397. set_j(1, 0, lmbd[0])
  398. set_j(1, 1, -lmbd[1])
  399. set_j(2, 1, lmbd[1])
  400. set_j(2, 2, -lmbd[2])
  401. return j
  402. def verify(self, zs, t):
  403. # Formulae derived by hand
  404. lmbd = np.array(self.lmbd)
  405. d10 = lmbd[1] - lmbd[0]
  406. d21 = lmbd[2] - lmbd[1]
  407. d20 = lmbd[2] - lmbd[0]
  408. e0 = np.exp(-lmbd[0] * t)
  409. e1 = np.exp(-lmbd[1] * t)
  410. e2 = np.exp(-lmbd[2] * t)
  411. u = np.vstack((
  412. self.z0[0] * e0,
  413. self.z0[1] * e1 + self.z0[0] * lmbd[0] / d10 * (e0 - e1),
  414. self.z0[2] * e2 + self.z0[1] * lmbd[1] / d21 * (e1 - e2) +
  415. lmbd[1] * lmbd[0] * self.z0[0] / d10 *
  416. (1 / d20 * (e0 - e2) - 1 / d21 * (e1 - e2)))).transpose()
  417. return allclose(u, zs, atol=self.atol, rtol=self.rtol)
  418. PROBLEMS = [SimpleOscillator, ComplexExp, Pi, CoupledDecay]
  419. #------------------------------------------------------------------------------
  420. def f(t, x):
  421. dxdt = [x[1], -x[0]]
  422. return dxdt
  423. def jac(t, x):
  424. j = array([[0.0, 1.0],
  425. [-1.0, 0.0]])
  426. return j
  427. def f1(t, x, omega):
  428. dxdt = [omega*x[1], -omega*x[0]]
  429. return dxdt
  430. def jac1(t, x, omega):
  431. j = array([[0.0, omega],
  432. [-omega, 0.0]])
  433. return j
  434. def f2(t, x, omega1, omega2):
  435. dxdt = [omega1*x[1], -omega2*x[0]]
  436. return dxdt
  437. def jac2(t, x, omega1, omega2):
  438. j = array([[0.0, omega1],
  439. [-omega2, 0.0]])
  440. return j
  441. def fv(t, x, omega):
  442. dxdt = [omega[0]*x[1], -omega[1]*x[0]]
  443. return dxdt
  444. def jacv(t, x, omega):
  445. j = array([[0.0, omega[0]],
  446. [-omega[1], 0.0]])
  447. return j
  448. class ODECheckParameterUse:
  449. """Call an ode-class solver with several cases of parameter use."""
  450. # solver_name must be set before tests can be run with this class.
  451. # Set these in subclasses.
  452. solver_name = ''
  453. solver_uses_jac = False
  454. def _get_solver(self, f, jac):
  455. solver = ode(f, jac)
  456. if self.solver_uses_jac:
  457. solver.set_integrator(self.solver_name, atol=1e-9, rtol=1e-7,
  458. with_jacobian=self.solver_uses_jac)
  459. else:
  460. # XXX Shouldn't set_integrator *always* accept the keyword arg
  461. # 'with_jacobian', and perhaps raise an exception if it is set
  462. # to True if the solver can't actually use it?
  463. solver.set_integrator(self.solver_name, atol=1e-9, rtol=1e-7)
  464. return solver
  465. def _check_solver(self, solver):
  466. ic = [1.0, 0.0]
  467. solver.set_initial_value(ic, 0.0)
  468. solver.integrate(pi)
  469. assert_array_almost_equal(solver.y, [-1.0, 0.0])
  470. def test_no_params(self):
  471. solver = self._get_solver(f, jac)
  472. self._check_solver(solver)
  473. def test_one_scalar_param(self):
  474. solver = self._get_solver(f1, jac1)
  475. omega = 1.0
  476. solver.set_f_params(omega)
  477. if self.solver_uses_jac:
  478. solver.set_jac_params(omega)
  479. self._check_solver(solver)
  480. def test_two_scalar_params(self):
  481. solver = self._get_solver(f2, jac2)
  482. omega1 = 1.0
  483. omega2 = 1.0
  484. solver.set_f_params(omega1, omega2)
  485. if self.solver_uses_jac:
  486. solver.set_jac_params(omega1, omega2)
  487. self._check_solver(solver)
  488. def test_vector_param(self):
  489. solver = self._get_solver(fv, jacv)
  490. omega = [1.0, 1.0]
  491. solver.set_f_params(omega)
  492. if self.solver_uses_jac:
  493. solver.set_jac_params(omega)
  494. self._check_solver(solver)
  495. def test_warns_on_failure(self):
  496. # Set nsteps small to ensure failure
  497. solver = self._get_solver(f, jac)
  498. solver.set_integrator(self.solver_name, nsteps=1)
  499. ic = [1.0, 0.0]
  500. solver.set_initial_value(ic, 0.0)
  501. assert_warns(UserWarning, solver.integrate, pi)
  502. class TestDOPRI5CheckParameterUse(ODECheckParameterUse):
  503. solver_name = 'dopri5'
  504. solver_uses_jac = False
  505. class TestDOP853CheckParameterUse(ODECheckParameterUse):
  506. solver_name = 'dop853'
  507. solver_uses_jac = False
  508. class TestVODECheckParameterUse(ODECheckParameterUse):
  509. solver_name = 'vode'
  510. solver_uses_jac = True
  511. class TestZVODECheckParameterUse(ODECheckParameterUse):
  512. solver_name = 'zvode'
  513. solver_uses_jac = True
  514. class TestLSODACheckParameterUse(ODECheckParameterUse):
  515. solver_name = 'lsoda'
  516. solver_uses_jac = True
  517. def test_odeint_trivial_time():
  518. # Test that odeint succeeds when given a single time point
  519. # and full_output=True. This is a regression test for gh-4282.
  520. y0 = 1
  521. t = [0]
  522. y, info = odeint(lambda y, t: -y, y0, t, full_output=True)
  523. assert_array_equal(y, np.array([[y0]]))
  524. def test_odeint_banded_jacobian():
  525. # Test the use of the `Dfun`, `ml` and `mu` options of odeint.
  526. def func(y, t, c):
  527. return c.dot(y)
  528. def jac(y, t, c):
  529. return c
  530. def jac_transpose(y, t, c):
  531. return c.T.copy(order='C')
  532. def bjac_rows(y, t, c):
  533. jac = np.row_stack((np.r_[0, np.diag(c, 1)],
  534. np.diag(c),
  535. np.r_[np.diag(c, -1), 0],
  536. np.r_[np.diag(c, -2), 0, 0]))
  537. return jac
  538. def bjac_cols(y, t, c):
  539. return bjac_rows(y, t, c).T.copy(order='C')
  540. c = array([[-205, 0.01, 0.00, 0.0],
  541. [0.1, -2.50, 0.02, 0.0],
  542. [1e-3, 0.01, -2.0, 0.01],
  543. [0.00, 0.00, 0.1, -1.0]])
  544. y0 = np.ones(4)
  545. t = np.array([0, 5, 10, 100])
  546. # Use the full Jacobian.
  547. sol1, info1 = odeint(func, y0, t, args=(c,), full_output=True,
  548. atol=1e-13, rtol=1e-11, mxstep=10000,
  549. Dfun=jac)
  550. # Use the transposed full Jacobian, with col_deriv=True.
  551. sol2, info2 = odeint(func, y0, t, args=(c,), full_output=True,
  552. atol=1e-13, rtol=1e-11, mxstep=10000,
  553. Dfun=jac_transpose, col_deriv=True)
  554. # Use the banded Jacobian.
  555. sol3, info3 = odeint(func, y0, t, args=(c,), full_output=True,
  556. atol=1e-13, rtol=1e-11, mxstep=10000,
  557. Dfun=bjac_rows, ml=2, mu=1)
  558. # Use the transposed banded Jacobian, with col_deriv=True.
  559. sol4, info4 = odeint(func, y0, t, args=(c,), full_output=True,
  560. atol=1e-13, rtol=1e-11, mxstep=10000,
  561. Dfun=bjac_cols, ml=2, mu=1, col_deriv=True)
  562. assert_allclose(sol1, sol2, err_msg="sol1 != sol2")
  563. assert_allclose(sol1, sol3, atol=1e-12, err_msg="sol1 != sol3")
  564. assert_allclose(sol3, sol4, err_msg="sol3 != sol4")
  565. # Verify that the number of jacobian evaluations was the same for the
  566. # calls of odeint with a full jacobian and with a banded jacobian. This is
  567. # a regression test--there was a bug in the handling of banded jacobians
  568. # that resulted in an incorrect jacobian matrix being passed to the LSODA
  569. # code. That would cause errors or excessive jacobian evaluations.
  570. assert_array_equal(info1['nje'], info2['nje'])
  571. assert_array_equal(info3['nje'], info4['nje'])
  572. # Test the use of tfirst
  573. sol1ty, info1ty = odeint(lambda t, y, c: func(y, t, c), y0, t, args=(c,),
  574. full_output=True, atol=1e-13, rtol=1e-11,
  575. mxstep=10000,
  576. Dfun=lambda t, y, c: jac(y, t, c), tfirst=True)
  577. # The code should execute the exact same sequence of floating point
  578. # calculations, so these should be exactly equal. We'll be safe and use
  579. # a small tolerance.
  580. assert_allclose(sol1, sol1ty, rtol=1e-12, err_msg="sol1 != sol1ty")
  581. def test_odeint_errors():
  582. def sys1d(x, t):
  583. return -100*x
  584. def bad1(x, t):
  585. return 1.0/0
  586. def bad2(x, t):
  587. return "foo"
  588. def bad_jac1(x, t):
  589. return 1.0/0
  590. def bad_jac2(x, t):
  591. return [["foo"]]
  592. def sys2d(x, t):
  593. return [-100*x[0], -0.1*x[1]]
  594. def sys2d_bad_jac(x, t):
  595. return [[1.0/0, 0], [0, -0.1]]
  596. assert_raises(ZeroDivisionError, odeint, bad1, 1.0, [0, 1])
  597. assert_raises(ValueError, odeint, bad2, 1.0, [0, 1])
  598. assert_raises(ZeroDivisionError, odeint, sys1d, 1.0, [0, 1], Dfun=bad_jac1)
  599. assert_raises(ValueError, odeint, sys1d, 1.0, [0, 1], Dfun=bad_jac2)
  600. assert_raises(ZeroDivisionError, odeint, sys2d, [1.0, 1.0], [0, 1],
  601. Dfun=sys2d_bad_jac)
  602. def test_odeint_bad_shapes():
  603. # Tests of some errors that can occur with odeint.
  604. def badrhs(x, t):
  605. return [1, -1]
  606. def sys1(x, t):
  607. return -100*x
  608. def badjac(x, t):
  609. return [[0, 0, 0]]
  610. # y0 must be at most 1-d.
  611. bad_y0 = [[0, 0], [0, 0]]
  612. assert_raises(ValueError, odeint, sys1, bad_y0, [0, 1])
  613. # t must be at most 1-d.
  614. bad_t = [[0, 1], [2, 3]]
  615. assert_raises(ValueError, odeint, sys1, [10.0], bad_t)
  616. # y0 is 10, but badrhs(x, t) returns [1, -1].
  617. assert_raises(RuntimeError, odeint, badrhs, 10, [0, 1])
  618. # shape of array returned by badjac(x, t) is not correct.
  619. assert_raises(RuntimeError, odeint, sys1, [10, 10], [0, 1], Dfun=badjac)
  620. def test_repeated_t_values():
  621. """Regression test for gh-8217."""
  622. def func(x, t):
  623. return -0.25*x
  624. t = np.zeros(10)
  625. sol = odeint(func, [1.], t)
  626. assert_array_equal(sol, np.ones((len(t), 1)))
  627. tau = 4*np.log(2)
  628. t = [0]*9 + [tau, 2*tau, 2*tau, 3*tau]
  629. sol = odeint(func, [1, 2], t, rtol=1e-12, atol=1e-12)
  630. expected_sol = np.array([[1.0, 2.0]]*9 +
  631. [[0.5, 1.0],
  632. [0.25, 0.5],
  633. [0.25, 0.5],
  634. [0.125, 0.25]])
  635. assert_allclose(sol, expected_sol)
  636. # Edge case: empty t sequence.
  637. sol = odeint(func, [1.], [])
  638. assert_array_equal(sol, np.array([], dtype=np.float64).reshape((0, 1)))
  639. # t values are not monotonic.
  640. assert_raises(ValueError, odeint, func, [1.], [0, 1, 0.5, 0])
  641. assert_raises(ValueError, odeint, func, [1, 2, 3], [0, -1, -2, 3])