test_bvp.py 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710
  1. import sys
  2. from io import StringIO
  3. import numpy as np
  4. from numpy.testing import (assert_, assert_array_equal, assert_allclose,
  5. assert_equal)
  6. from pytest import raises as assert_raises
  7. from scipy.sparse import coo_matrix
  8. from scipy.special import erf
  9. from scipy.integrate._bvp import (modify_mesh, estimate_fun_jac,
  10. estimate_bc_jac, compute_jac_indices,
  11. construct_global_jac, solve_bvp)
  12. import pytest
  13. def exp_fun(x, y):
  14. return np.vstack((y[1], y[0]))
  15. def exp_fun_jac(x, y):
  16. df_dy = np.empty((2, 2, x.shape[0]))
  17. df_dy[0, 0] = 0
  18. df_dy[0, 1] = 1
  19. df_dy[1, 0] = 1
  20. df_dy[1, 1] = 0
  21. return df_dy
  22. def exp_bc(ya, yb):
  23. return np.hstack((ya[0] - 1, yb[0]))
  24. def exp_bc_complex(ya, yb):
  25. return np.hstack((ya[0] - 1 - 1j, yb[0]))
  26. def exp_bc_jac(ya, yb):
  27. dbc_dya = np.array([
  28. [1, 0],
  29. [0, 0]
  30. ])
  31. dbc_dyb = np.array([
  32. [0, 0],
  33. [1, 0]
  34. ])
  35. return dbc_dya, dbc_dyb
  36. def exp_sol(x):
  37. return (np.exp(-x) - np.exp(x - 2)) / (1 - np.exp(-2))
  38. def sl_fun(x, y, p):
  39. return np.vstack((y[1], -p[0]**2 * y[0]))
  40. def sl_fun_jac(x, y, p):
  41. n, m = y.shape
  42. df_dy = np.empty((n, 2, m))
  43. df_dy[0, 0] = 0
  44. df_dy[0, 1] = 1
  45. df_dy[1, 0] = -p[0]**2
  46. df_dy[1, 1] = 0
  47. df_dp = np.empty((n, 1, m))
  48. df_dp[0, 0] = 0
  49. df_dp[1, 0] = -2 * p[0] * y[0]
  50. return df_dy, df_dp
  51. def sl_bc(ya, yb, p):
  52. return np.hstack((ya[0], yb[0], ya[1] - p[0]))
  53. def sl_bc_jac(ya, yb, p):
  54. dbc_dya = np.zeros((3, 2))
  55. dbc_dya[0, 0] = 1
  56. dbc_dya[2, 1] = 1
  57. dbc_dyb = np.zeros((3, 2))
  58. dbc_dyb[1, 0] = 1
  59. dbc_dp = np.zeros((3, 1))
  60. dbc_dp[2, 0] = -1
  61. return dbc_dya, dbc_dyb, dbc_dp
  62. def sl_sol(x, p):
  63. return np.sin(p[0] * x)
  64. def emden_fun(x, y):
  65. return np.vstack((y[1], -y[0]**5))
  66. def emden_fun_jac(x, y):
  67. df_dy = np.empty((2, 2, x.shape[0]))
  68. df_dy[0, 0] = 0
  69. df_dy[0, 1] = 1
  70. df_dy[1, 0] = -5 * y[0]**4
  71. df_dy[1, 1] = 0
  72. return df_dy
  73. def emden_bc(ya, yb):
  74. return np.array([ya[1], yb[0] - (3/4)**0.5])
  75. def emden_bc_jac(ya, yb):
  76. dbc_dya = np.array([
  77. [0, 1],
  78. [0, 0]
  79. ])
  80. dbc_dyb = np.array([
  81. [0, 0],
  82. [1, 0]
  83. ])
  84. return dbc_dya, dbc_dyb
  85. def emden_sol(x):
  86. return (1 + x**2/3)**-0.5
  87. def undefined_fun(x, y):
  88. return np.zeros_like(y)
  89. def undefined_bc(ya, yb):
  90. return np.array([ya[0], yb[0] - 1])
  91. def big_fun(x, y):
  92. f = np.zeros_like(y)
  93. f[::2] = y[1::2]
  94. return f
  95. def big_bc(ya, yb):
  96. return np.hstack((ya[::2], yb[::2] - 1))
  97. def big_sol(x, n):
  98. y = np.ones((2 * n, x.size))
  99. y[::2] = x
  100. return x
  101. def big_fun_with_parameters(x, y, p):
  102. """ Big version of sl_fun, with two parameters.
  103. The two differential equations represented by sl_fun are broadcast to the
  104. number of rows of y, rotating between the parameters p[0] and p[1].
  105. Here are the differential equations:
  106. dy[0]/dt = y[1]
  107. dy[1]/dt = -p[0]**2 * y[0]
  108. dy[2]/dt = y[3]
  109. dy[3]/dt = -p[1]**2 * y[2]
  110. dy[4]/dt = y[5]
  111. dy[5]/dt = -p[0]**2 * y[4]
  112. dy[6]/dt = y[7]
  113. dy[7]/dt = -p[1]**2 * y[6]
  114. .
  115. .
  116. .
  117. """
  118. f = np.zeros_like(y)
  119. f[::2] = y[1::2]
  120. f[1::4] = -p[0]**2 * y[::4]
  121. f[3::4] = -p[1]**2 * y[2::4]
  122. return f
  123. def big_fun_with_parameters_jac(x, y, p):
  124. # big version of sl_fun_jac, with two parameters
  125. n, m = y.shape
  126. df_dy = np.zeros((n, n, m))
  127. df_dy[range(0, n, 2), range(1, n, 2)] = 1
  128. df_dy[range(1, n, 4), range(0, n, 4)] = -p[0]**2
  129. df_dy[range(3, n, 4), range(2, n, 4)] = -p[1]**2
  130. df_dp = np.zeros((n, 2, m))
  131. df_dp[range(1, n, 4), 0] = -2 * p[0] * y[range(0, n, 4)]
  132. df_dp[range(3, n, 4), 1] = -2 * p[1] * y[range(2, n, 4)]
  133. return df_dy, df_dp
  134. def big_bc_with_parameters(ya, yb, p):
  135. # big version of sl_bc, with two parameters
  136. return np.hstack((ya[::2], yb[::2], ya[1] - p[0], ya[3] - p[1]))
  137. def big_bc_with_parameters_jac(ya, yb, p):
  138. # big version of sl_bc_jac, with two parameters
  139. n = ya.shape[0]
  140. dbc_dya = np.zeros((n + 2, n))
  141. dbc_dyb = np.zeros((n + 2, n))
  142. dbc_dya[range(n // 2), range(0, n, 2)] = 1
  143. dbc_dyb[range(n // 2, n), range(0, n, 2)] = 1
  144. dbc_dp = np.zeros((n + 2, 2))
  145. dbc_dp[n, 0] = -1
  146. dbc_dya[n, 1] = 1
  147. dbc_dp[n + 1, 1] = -1
  148. dbc_dya[n + 1, 3] = 1
  149. return dbc_dya, dbc_dyb, dbc_dp
  150. def big_sol_with_parameters(x, p):
  151. # big version of sl_sol, with two parameters
  152. return np.vstack((np.sin(p[0] * x), np.sin(p[1] * x)))
  153. def shock_fun(x, y):
  154. eps = 1e-3
  155. return np.vstack((
  156. y[1],
  157. -(x * y[1] + eps * np.pi**2 * np.cos(np.pi * x) +
  158. np.pi * x * np.sin(np.pi * x)) / eps
  159. ))
  160. def shock_bc(ya, yb):
  161. return np.array([ya[0] + 2, yb[0]])
  162. def shock_sol(x):
  163. eps = 1e-3
  164. k = np.sqrt(2 * eps)
  165. return np.cos(np.pi * x) + erf(x / k) / erf(1 / k)
  166. def nonlin_bc_fun(x, y):
  167. # laplace eq.
  168. return np.stack([y[1], np.zeros_like(x)])
  169. def nonlin_bc_bc(ya, yb):
  170. phiA, phipA = ya
  171. phiC, phipC = yb
  172. kappa, ioA, ioC, V, f = 1.64, 0.01, 1.0e-4, 0.5, 38.9
  173. # Butler-Volmer Kinetics at Anode
  174. hA = 0.0-phiA-0.0
  175. iA = ioA * (np.exp(f*hA) - np.exp(-f*hA))
  176. res0 = iA + kappa * phipA
  177. # Butler-Volmer Kinetics at Cathode
  178. hC = V - phiC - 1.0
  179. iC = ioC * (np.exp(f*hC) - np.exp(-f*hC))
  180. res1 = iC - kappa*phipC
  181. return np.array([res0, res1])
  182. def nonlin_bc_sol(x):
  183. return -0.13426436116763119 - 1.1308709 * x
  184. def test_modify_mesh():
  185. x = np.array([0, 1, 3, 9], dtype=float)
  186. x_new = modify_mesh(x, np.array([0]), np.array([2]))
  187. assert_array_equal(x_new, np.array([0, 0.5, 1, 3, 5, 7, 9]))
  188. x = np.array([-6, -3, 0, 3, 6], dtype=float)
  189. x_new = modify_mesh(x, np.array([1], dtype=int), np.array([0, 2, 3]))
  190. assert_array_equal(x_new, [-6, -5, -4, -3, -1.5, 0, 1, 2, 3, 4, 5, 6])
  191. def test_compute_fun_jac():
  192. x = np.linspace(0, 1, 5)
  193. y = np.empty((2, x.shape[0]))
  194. y[0] = 0.01
  195. y[1] = 0.02
  196. p = np.array([])
  197. df_dy, df_dp = estimate_fun_jac(lambda x, y, p: exp_fun(x, y), x, y, p)
  198. df_dy_an = exp_fun_jac(x, y)
  199. assert_allclose(df_dy, df_dy_an)
  200. assert_(df_dp is None)
  201. x = np.linspace(0, np.pi, 5)
  202. y = np.empty((2, x.shape[0]))
  203. y[0] = np.sin(x)
  204. y[1] = np.cos(x)
  205. p = np.array([1.0])
  206. df_dy, df_dp = estimate_fun_jac(sl_fun, x, y, p)
  207. df_dy_an, df_dp_an = sl_fun_jac(x, y, p)
  208. assert_allclose(df_dy, df_dy_an)
  209. assert_allclose(df_dp, df_dp_an)
  210. x = np.linspace(0, 1, 10)
  211. y = np.empty((2, x.shape[0]))
  212. y[0] = (3/4)**0.5
  213. y[1] = 1e-4
  214. p = np.array([])
  215. df_dy, df_dp = estimate_fun_jac(lambda x, y, p: emden_fun(x, y), x, y, p)
  216. df_dy_an = emden_fun_jac(x, y)
  217. assert_allclose(df_dy, df_dy_an)
  218. assert_(df_dp is None)
  219. def test_compute_bc_jac():
  220. ya = np.array([-1.0, 2])
  221. yb = np.array([0.5, 3])
  222. p = np.array([])
  223. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(
  224. lambda ya, yb, p: exp_bc(ya, yb), ya, yb, p)
  225. dbc_dya_an, dbc_dyb_an = exp_bc_jac(ya, yb)
  226. assert_allclose(dbc_dya, dbc_dya_an)
  227. assert_allclose(dbc_dyb, dbc_dyb_an)
  228. assert_(dbc_dp is None)
  229. ya = np.array([0.0, 1])
  230. yb = np.array([0.0, -1])
  231. p = np.array([0.5])
  232. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(sl_bc, ya, yb, p)
  233. dbc_dya_an, dbc_dyb_an, dbc_dp_an = sl_bc_jac(ya, yb, p)
  234. assert_allclose(dbc_dya, dbc_dya_an)
  235. assert_allclose(dbc_dyb, dbc_dyb_an)
  236. assert_allclose(dbc_dp, dbc_dp_an)
  237. ya = np.array([0.5, 100])
  238. yb = np.array([-1000, 10.5])
  239. p = np.array([])
  240. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(
  241. lambda ya, yb, p: emden_bc(ya, yb), ya, yb, p)
  242. dbc_dya_an, dbc_dyb_an = emden_bc_jac(ya, yb)
  243. assert_allclose(dbc_dya, dbc_dya_an)
  244. assert_allclose(dbc_dyb, dbc_dyb_an)
  245. assert_(dbc_dp is None)
  246. def test_compute_jac_indices():
  247. n = 2
  248. m = 4
  249. k = 2
  250. i, j = compute_jac_indices(n, m, k)
  251. s = coo_matrix((np.ones_like(i), (i, j))).toarray()
  252. s_true = np.array([
  253. [1, 1, 1, 1, 0, 0, 0, 0, 1, 1],
  254. [1, 1, 1, 1, 0, 0, 0, 0, 1, 1],
  255. [0, 0, 1, 1, 1, 1, 0, 0, 1, 1],
  256. [0, 0, 1, 1, 1, 1, 0, 0, 1, 1],
  257. [0, 0, 0, 0, 1, 1, 1, 1, 1, 1],
  258. [0, 0, 0, 0, 1, 1, 1, 1, 1, 1],
  259. [1, 1, 0, 0, 0, 0, 1, 1, 1, 1],
  260. [1, 1, 0, 0, 0, 0, 1, 1, 1, 1],
  261. [1, 1, 0, 0, 0, 0, 1, 1, 1, 1],
  262. [1, 1, 0, 0, 0, 0, 1, 1, 1, 1],
  263. ])
  264. assert_array_equal(s, s_true)
  265. def test_compute_global_jac():
  266. n = 2
  267. m = 5
  268. k = 1
  269. i_jac, j_jac = compute_jac_indices(2, 5, 1)
  270. x = np.linspace(0, 1, 5)
  271. h = np.diff(x)
  272. y = np.vstack((np.sin(np.pi * x), np.pi * np.cos(np.pi * x)))
  273. p = np.array([3.0])
  274. f = sl_fun(x, y, p)
  275. x_middle = x[:-1] + 0.5 * h
  276. y_middle = 0.5 * (y[:, :-1] + y[:, 1:]) - h/8 * (f[:, 1:] - f[:, :-1])
  277. df_dy, df_dp = sl_fun_jac(x, y, p)
  278. df_dy_middle, df_dp_middle = sl_fun_jac(x_middle, y_middle, p)
  279. dbc_dya, dbc_dyb, dbc_dp = sl_bc_jac(y[:, 0], y[:, -1], p)
  280. J = construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy, df_dy_middle,
  281. df_dp, df_dp_middle, dbc_dya, dbc_dyb, dbc_dp)
  282. J = J.toarray()
  283. def J_block(h, p):
  284. return np.array([
  285. [h**2*p**2/12 - 1, -0.5*h, -h**2*p**2/12 + 1, -0.5*h],
  286. [0.5*h*p**2, h**2*p**2/12 - 1, 0.5*h*p**2, 1 - h**2*p**2/12]
  287. ])
  288. J_true = np.zeros((m * n + k, m * n + k))
  289. for i in range(m - 1):
  290. J_true[i * n: (i + 1) * n, i * n: (i + 2) * n] = J_block(h[i], p[0])
  291. J_true[:(m - 1) * n:2, -1] = p * h**2/6 * (y[0, :-1] - y[0, 1:])
  292. J_true[1:(m - 1) * n:2, -1] = p * (h * (y[0, :-1] + y[0, 1:]) +
  293. h**2/6 * (y[1, :-1] - y[1, 1:]))
  294. J_true[8, 0] = 1
  295. J_true[9, 8] = 1
  296. J_true[10, 1] = 1
  297. J_true[10, 10] = -1
  298. assert_allclose(J, J_true, rtol=1e-10)
  299. df_dy, df_dp = estimate_fun_jac(sl_fun, x, y, p)
  300. df_dy_middle, df_dp_middle = estimate_fun_jac(sl_fun, x_middle, y_middle, p)
  301. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(sl_bc, y[:, 0], y[:, -1], p)
  302. J = construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy, df_dy_middle,
  303. df_dp, df_dp_middle, dbc_dya, dbc_dyb, dbc_dp)
  304. J = J.toarray()
  305. assert_allclose(J, J_true, rtol=2e-8, atol=2e-8)
  306. def test_parameter_validation():
  307. x = [0, 1, 0.5]
  308. y = np.zeros((2, 3))
  309. assert_raises(ValueError, solve_bvp, exp_fun, exp_bc, x, y)
  310. x = np.linspace(0, 1, 5)
  311. y = np.zeros((2, 4))
  312. assert_raises(ValueError, solve_bvp, exp_fun, exp_bc, x, y)
  313. def fun(x, y, p):
  314. return exp_fun(x, y)
  315. def bc(ya, yb, p):
  316. return exp_bc(ya, yb)
  317. y = np.zeros((2, x.shape[0]))
  318. assert_raises(ValueError, solve_bvp, fun, bc, x, y, p=[1])
  319. def wrong_shape_fun(x, y):
  320. return np.zeros(3)
  321. assert_raises(ValueError, solve_bvp, wrong_shape_fun, bc, x, y)
  322. S = np.array([[0, 0]])
  323. assert_raises(ValueError, solve_bvp, exp_fun, exp_bc, x, y, S=S)
  324. def test_no_params():
  325. x = np.linspace(0, 1, 5)
  326. x_test = np.linspace(0, 1, 100)
  327. y = np.zeros((2, x.shape[0]))
  328. for fun_jac in [None, exp_fun_jac]:
  329. for bc_jac in [None, exp_bc_jac]:
  330. sol = solve_bvp(exp_fun, exp_bc, x, y, fun_jac=fun_jac,
  331. bc_jac=bc_jac)
  332. assert_equal(sol.status, 0)
  333. assert_(sol.success)
  334. assert_equal(sol.x.size, 5)
  335. sol_test = sol.sol(x_test)
  336. assert_allclose(sol_test[0], exp_sol(x_test), atol=1e-5)
  337. f_test = exp_fun(x_test, sol_test)
  338. r = sol.sol(x_test, 1) - f_test
  339. rel_res = r / (1 + np.abs(f_test))
  340. norm_res = np.sum(rel_res**2, axis=0)**0.5
  341. assert_(np.all(norm_res < 1e-3))
  342. assert_(np.all(sol.rms_residuals < 1e-3))
  343. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  344. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  345. def test_with_params():
  346. x = np.linspace(0, np.pi, 5)
  347. x_test = np.linspace(0, np.pi, 100)
  348. y = np.ones((2, x.shape[0]))
  349. for fun_jac in [None, sl_fun_jac]:
  350. for bc_jac in [None, sl_bc_jac]:
  351. sol = solve_bvp(sl_fun, sl_bc, x, y, p=[0.5], fun_jac=fun_jac,
  352. bc_jac=bc_jac)
  353. assert_equal(sol.status, 0)
  354. assert_(sol.success)
  355. assert_(sol.x.size < 10)
  356. assert_allclose(sol.p, [1], rtol=1e-4)
  357. sol_test = sol.sol(x_test)
  358. assert_allclose(sol_test[0], sl_sol(x_test, [1]),
  359. rtol=1e-4, atol=1e-4)
  360. f_test = sl_fun(x_test, sol_test, [1])
  361. r = sol.sol(x_test, 1) - f_test
  362. rel_res = r / (1 + np.abs(f_test))
  363. norm_res = np.sum(rel_res ** 2, axis=0) ** 0.5
  364. assert_(np.all(norm_res < 1e-3))
  365. assert_(np.all(sol.rms_residuals < 1e-3))
  366. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  367. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  368. def test_singular_term():
  369. x = np.linspace(0, 1, 10)
  370. x_test = np.linspace(0.05, 1, 100)
  371. y = np.empty((2, 10))
  372. y[0] = (3/4)**0.5
  373. y[1] = 1e-4
  374. S = np.array([[0, 0], [0, -2]])
  375. for fun_jac in [None, emden_fun_jac]:
  376. for bc_jac in [None, emden_bc_jac]:
  377. sol = solve_bvp(emden_fun, emden_bc, x, y, S=S, fun_jac=fun_jac,
  378. bc_jac=bc_jac)
  379. assert_equal(sol.status, 0)
  380. assert_(sol.success)
  381. assert_equal(sol.x.size, 10)
  382. sol_test = sol.sol(x_test)
  383. assert_allclose(sol_test[0], emden_sol(x_test), atol=1e-5)
  384. f_test = emden_fun(x_test, sol_test) + S.dot(sol_test) / x_test
  385. r = sol.sol(x_test, 1) - f_test
  386. rel_res = r / (1 + np.abs(f_test))
  387. norm_res = np.sum(rel_res ** 2, axis=0) ** 0.5
  388. assert_(np.all(norm_res < 1e-3))
  389. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  390. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  391. def test_complex():
  392. # The test is essentially the same as test_no_params, but boundary
  393. # conditions are turned into complex.
  394. x = np.linspace(0, 1, 5)
  395. x_test = np.linspace(0, 1, 100)
  396. y = np.zeros((2, x.shape[0]), dtype=complex)
  397. for fun_jac in [None, exp_fun_jac]:
  398. for bc_jac in [None, exp_bc_jac]:
  399. sol = solve_bvp(exp_fun, exp_bc_complex, x, y, fun_jac=fun_jac,
  400. bc_jac=bc_jac)
  401. assert_equal(sol.status, 0)
  402. assert_(sol.success)
  403. sol_test = sol.sol(x_test)
  404. assert_allclose(sol_test[0].real, exp_sol(x_test), atol=1e-5)
  405. assert_allclose(sol_test[0].imag, exp_sol(x_test), atol=1e-5)
  406. f_test = exp_fun(x_test, sol_test)
  407. r = sol.sol(x_test, 1) - f_test
  408. rel_res = r / (1 + np.abs(f_test))
  409. norm_res = np.sum(np.real(rel_res * np.conj(rel_res)),
  410. axis=0) ** 0.5
  411. assert_(np.all(norm_res < 1e-3))
  412. assert_(np.all(sol.rms_residuals < 1e-3))
  413. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  414. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  415. def test_failures():
  416. x = np.linspace(0, 1, 2)
  417. y = np.zeros((2, x.size))
  418. res = solve_bvp(exp_fun, exp_bc, x, y, tol=1e-5, max_nodes=5)
  419. assert_equal(res.status, 1)
  420. assert_(not res.success)
  421. x = np.linspace(0, 1, 5)
  422. y = np.zeros((2, x.size))
  423. res = solve_bvp(undefined_fun, undefined_bc, x, y)
  424. assert_equal(res.status, 2)
  425. assert_(not res.success)
  426. def test_big_problem():
  427. n = 30
  428. x = np.linspace(0, 1, 5)
  429. y = np.zeros((2 * n, x.size))
  430. sol = solve_bvp(big_fun, big_bc, x, y)
  431. assert_equal(sol.status, 0)
  432. assert_(sol.success)
  433. sol_test = sol.sol(x)
  434. assert_allclose(sol_test[0], big_sol(x, n))
  435. f_test = big_fun(x, sol_test)
  436. r = sol.sol(x, 1) - f_test
  437. rel_res = r / (1 + np.abs(f_test))
  438. norm_res = np.sum(np.real(rel_res * np.conj(rel_res)), axis=0) ** 0.5
  439. assert_(np.all(norm_res < 1e-3))
  440. assert_(np.all(sol.rms_residuals < 1e-3))
  441. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  442. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  443. def test_big_problem_with_parameters():
  444. n = 30
  445. x = np.linspace(0, np.pi, 5)
  446. x_test = np.linspace(0, np.pi, 100)
  447. y = np.ones((2 * n, x.size))
  448. for fun_jac in [None, big_fun_with_parameters_jac]:
  449. for bc_jac in [None, big_bc_with_parameters_jac]:
  450. sol = solve_bvp(big_fun_with_parameters, big_bc_with_parameters, x,
  451. y, p=[0.5, 0.5], fun_jac=fun_jac, bc_jac=bc_jac)
  452. assert_equal(sol.status, 0)
  453. assert_(sol.success)
  454. assert_allclose(sol.p, [1, 1], rtol=1e-4)
  455. sol_test = sol.sol(x_test)
  456. for isol in range(0, n, 4):
  457. assert_allclose(sol_test[isol],
  458. big_sol_with_parameters(x_test, [1, 1])[0],
  459. rtol=1e-4, atol=1e-4)
  460. assert_allclose(sol_test[isol + 2],
  461. big_sol_with_parameters(x_test, [1, 1])[1],
  462. rtol=1e-4, atol=1e-4)
  463. f_test = big_fun_with_parameters(x_test, sol_test, [1, 1])
  464. r = sol.sol(x_test, 1) - f_test
  465. rel_res = r / (1 + np.abs(f_test))
  466. norm_res = np.sum(rel_res ** 2, axis=0) ** 0.5
  467. assert_(np.all(norm_res < 1e-3))
  468. assert_(np.all(sol.rms_residuals < 1e-3))
  469. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  470. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  471. def test_shock_layer():
  472. x = np.linspace(-1, 1, 5)
  473. x_test = np.linspace(-1, 1, 100)
  474. y = np.zeros((2, x.size))
  475. sol = solve_bvp(shock_fun, shock_bc, x, y)
  476. assert_equal(sol.status, 0)
  477. assert_(sol.success)
  478. assert_(sol.x.size < 110)
  479. sol_test = sol.sol(x_test)
  480. assert_allclose(sol_test[0], shock_sol(x_test), rtol=1e-5, atol=1e-5)
  481. f_test = shock_fun(x_test, sol_test)
  482. r = sol.sol(x_test, 1) - f_test
  483. rel_res = r / (1 + np.abs(f_test))
  484. norm_res = np.sum(rel_res ** 2, axis=0) ** 0.5
  485. assert_(np.all(norm_res < 1e-3))
  486. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  487. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  488. def test_nonlin_bc():
  489. x = np.linspace(0, 0.1, 5)
  490. x_test = x
  491. y = np.zeros([2, x.size])
  492. sol = solve_bvp(nonlin_bc_fun, nonlin_bc_bc, x, y)
  493. assert_equal(sol.status, 0)
  494. assert_(sol.success)
  495. assert_(sol.x.size < 8)
  496. sol_test = sol.sol(x_test)
  497. assert_allclose(sol_test[0], nonlin_bc_sol(x_test), rtol=1e-5, atol=1e-5)
  498. f_test = nonlin_bc_fun(x_test, sol_test)
  499. r = sol.sol(x_test, 1) - f_test
  500. rel_res = r / (1 + np.abs(f_test))
  501. norm_res = np.sum(rel_res ** 2, axis=0) ** 0.5
  502. assert_(np.all(norm_res < 1e-3))
  503. assert_allclose(sol.sol(sol.x), sol.y, rtol=1e-10, atol=1e-10)
  504. assert_allclose(sol.sol(sol.x, 1), sol.yp, rtol=1e-10, atol=1e-10)
  505. @pytest.mark.thread_unsafe(reason="multithreaded sys.stdout parsing is not thread-safe")
  506. def test_verbose():
  507. # Smoke test that checks the printing does something and does not crash
  508. x = np.linspace(0, 1, 5)
  509. y = np.zeros((2, x.shape[0]))
  510. for verbose in [0, 1, 2]:
  511. old_stdout = sys.stdout
  512. sys.stdout = StringIO()
  513. try:
  514. sol = solve_bvp(exp_fun, exp_bc, x, y, verbose=verbose)
  515. text = sys.stdout.getvalue()
  516. finally:
  517. sys.stdout = old_stdout
  518. assert_(sol.success)
  519. if verbose == 0:
  520. assert_(not text, text)
  521. if verbose >= 1:
  522. assert_("Solved in" in text, text)
  523. if verbose >= 2:
  524. assert_("Max residual" in text, text)