_bvp.py 40 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162
  1. """Boundary value problem solver."""
  2. from warnings import warn
  3. import numpy as np
  4. from numpy.linalg import pinv
  5. from scipy.sparse import coo_matrix, csc_matrix
  6. from scipy.sparse.linalg import splu
  7. from scipy.optimize import OptimizeResult
  8. from scipy._lib._array_api import xp_capabilities
  9. EPS = np.finfo(float).eps
  10. def estimate_fun_jac(fun, x, y, p, f0=None):
  11. """Estimate derivatives of an ODE system rhs with forward differences.
  12. Returns
  13. -------
  14. df_dy : ndarray, shape (n, n, m)
  15. Derivatives with respect to y. An element (i, j, q) corresponds to
  16. d f_i(x_q, y_q) / d (y_q)_j.
  17. df_dp : ndarray with shape (n, k, m) or None
  18. Derivatives with respect to p. An element (i, j, q) corresponds to
  19. d f_i(x_q, y_q, p) / d p_j. If `p` is empty, None is returned.
  20. """
  21. n, m = y.shape
  22. if f0 is None:
  23. f0 = fun(x, y, p)
  24. dtype = y.dtype
  25. df_dy = np.empty((n, n, m), dtype=dtype)
  26. h = EPS**0.5 * (1 + np.abs(y))
  27. for i in range(n):
  28. y_new = y.copy()
  29. y_new[i] += h[i]
  30. hi = y_new[i] - y[i]
  31. f_new = fun(x, y_new, p)
  32. df_dy[:, i, :] = (f_new - f0) / hi
  33. k = p.shape[0]
  34. if k == 0:
  35. df_dp = None
  36. else:
  37. df_dp = np.empty((n, k, m), dtype=dtype)
  38. h = EPS**0.5 * (1 + np.abs(p))
  39. for i in range(k):
  40. p_new = p.copy()
  41. p_new[i] += h[i]
  42. hi = p_new[i] - p[i]
  43. f_new = fun(x, y, p_new)
  44. df_dp[:, i, :] = (f_new - f0) / hi
  45. return df_dy, df_dp
  46. def estimate_bc_jac(bc, ya, yb, p, bc0=None):
  47. """Estimate derivatives of boundary conditions with forward differences.
  48. Returns
  49. -------
  50. dbc_dya : ndarray, shape (n + k, n)
  51. Derivatives with respect to ya. An element (i, j) corresponds to
  52. d bc_i / d ya_j.
  53. dbc_dyb : ndarray, shape (n + k, n)
  54. Derivatives with respect to yb. An element (i, j) corresponds to
  55. d bc_i / d ya_j.
  56. dbc_dp : ndarray with shape (n + k, k) or None
  57. Derivatives with respect to p. An element (i, j) corresponds to
  58. d bc_i / d p_j. If `p` is empty, None is returned.
  59. """
  60. n = ya.shape[0]
  61. k = p.shape[0]
  62. if bc0 is None:
  63. bc0 = bc(ya, yb, p)
  64. dtype = ya.dtype
  65. dbc_dya = np.empty((n, n + k), dtype=dtype)
  66. h = EPS**0.5 * (1 + np.abs(ya))
  67. for i in range(n):
  68. ya_new = ya.copy()
  69. ya_new[i] += h[i]
  70. hi = ya_new[i] - ya[i]
  71. bc_new = bc(ya_new, yb, p)
  72. dbc_dya[i] = (bc_new - bc0) / hi
  73. dbc_dya = dbc_dya.T
  74. h = EPS**0.5 * (1 + np.abs(yb))
  75. dbc_dyb = np.empty((n, n + k), dtype=dtype)
  76. for i in range(n):
  77. yb_new = yb.copy()
  78. yb_new[i] += h[i]
  79. hi = yb_new[i] - yb[i]
  80. bc_new = bc(ya, yb_new, p)
  81. dbc_dyb[i] = (bc_new - bc0) / hi
  82. dbc_dyb = dbc_dyb.T
  83. if k == 0:
  84. dbc_dp = None
  85. else:
  86. h = EPS**0.5 * (1 + np.abs(p))
  87. dbc_dp = np.empty((k, n + k), dtype=dtype)
  88. for i in range(k):
  89. p_new = p.copy()
  90. p_new[i] += h[i]
  91. hi = p_new[i] - p[i]
  92. bc_new = bc(ya, yb, p_new)
  93. dbc_dp[i] = (bc_new - bc0) / hi
  94. dbc_dp = dbc_dp.T
  95. return dbc_dya, dbc_dyb, dbc_dp
  96. def compute_jac_indices(n, m, k):
  97. """Compute indices for the collocation system Jacobian construction.
  98. See `construct_global_jac` for the explanation.
  99. """
  100. i_col = np.repeat(np.arange((m - 1) * n), n)
  101. j_col = (np.tile(np.arange(n), n * (m - 1)) +
  102. np.repeat(np.arange(m - 1) * n, n**2))
  103. i_bc = np.repeat(np.arange((m - 1) * n, m * n + k), n)
  104. j_bc = np.tile(np.arange(n), n + k)
  105. i_p_col = np.repeat(np.arange((m - 1) * n), k)
  106. j_p_col = np.tile(np.arange(m * n, m * n + k), (m - 1) * n)
  107. i_p_bc = np.repeat(np.arange((m - 1) * n, m * n + k), k)
  108. j_p_bc = np.tile(np.arange(m * n, m * n + k), n + k)
  109. i = np.hstack((i_col, i_col, i_bc, i_bc, i_p_col, i_p_bc))
  110. j = np.hstack((j_col, j_col + n,
  111. j_bc, j_bc + (m - 1) * n,
  112. j_p_col, j_p_bc))
  113. return i, j
  114. def stacked_matmul(a, b):
  115. """Stacked matrix multiply: out[i,:,:] = np.dot(a[i,:,:], b[i,:,:]).
  116. Empirical optimization. Use outer Python loop and BLAS for large
  117. matrices, otherwise use a single einsum call.
  118. """
  119. if a.shape[1] > 50:
  120. out = np.empty((a.shape[0], a.shape[1], b.shape[2]))
  121. for i in range(a.shape[0]):
  122. out[i] = np.dot(a[i], b[i])
  123. return out
  124. else:
  125. return np.einsum('...ij,...jk->...ik', a, b)
  126. def construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy, df_dy_middle, df_dp,
  127. df_dp_middle, dbc_dya, dbc_dyb, dbc_dp):
  128. """Construct the Jacobian of the collocation system.
  129. There are n * m + k functions: m - 1 collocations residuals, each
  130. containing n components, followed by n + k boundary condition residuals.
  131. There are n * m + k variables: m vectors of y, each containing n
  132. components, followed by k values of vector p.
  133. For example, let m = 4, n = 2 and k = 1, then the Jacobian will have
  134. the following sparsity structure:
  135. 1 1 2 2 0 0 0 0 5
  136. 1 1 2 2 0 0 0 0 5
  137. 0 0 1 1 2 2 0 0 5
  138. 0 0 1 1 2 2 0 0 5
  139. 0 0 0 0 1 1 2 2 5
  140. 0 0 0 0 1 1 2 2 5
  141. 3 3 0 0 0 0 4 4 6
  142. 3 3 0 0 0 0 4 4 6
  143. 3 3 0 0 0 0 4 4 6
  144. Zeros denote identically zero values, other values denote different kinds
  145. of blocks in the matrix (see below). The blank row indicates the separation
  146. of collocation residuals from boundary conditions. And the blank column
  147. indicates the separation of y values from p values.
  148. Refer to [1]_ (p. 306) for the formula of n x n blocks for derivatives
  149. of collocation residuals with respect to y.
  150. Parameters
  151. ----------
  152. n : int
  153. Number of equations in the ODE system.
  154. m : int
  155. Number of nodes in the mesh.
  156. k : int
  157. Number of the unknown parameters.
  158. i_jac, j_jac : ndarray
  159. Row and column indices returned by `compute_jac_indices`. They
  160. represent different blocks in the Jacobian matrix in the following
  161. order (see the scheme above):
  162. * 1: m - 1 diagonal n x n blocks for the collocation residuals.
  163. * 2: m - 1 off-diagonal n x n blocks for the collocation residuals.
  164. * 3 : (n + k) x n block for the dependency of the boundary
  165. conditions on ya.
  166. * 4: (n + k) x n block for the dependency of the boundary
  167. conditions on yb.
  168. * 5: (m - 1) * n x k block for the dependency of the collocation
  169. residuals on p.
  170. * 6: (n + k) x k block for the dependency of the boundary
  171. conditions on p.
  172. df_dy : ndarray, shape (n, n, m)
  173. Jacobian of f with respect to y computed at the mesh nodes.
  174. df_dy_middle : ndarray, shape (n, n, m - 1)
  175. Jacobian of f with respect to y computed at the middle between the
  176. mesh nodes.
  177. df_dp : ndarray with shape (n, k, m) or None
  178. Jacobian of f with respect to p computed at the mesh nodes.
  179. df_dp_middle : ndarray with shape (n, k, m - 1) or None
  180. Jacobian of f with respect to p computed at the middle between the
  181. mesh nodes.
  182. dbc_dya, dbc_dyb : ndarray, shape (n, n)
  183. Jacobian of bc with respect to ya and yb.
  184. dbc_dp : ndarray with shape (n, k) or None
  185. Jacobian of bc with respect to p.
  186. Returns
  187. -------
  188. J : csc_matrix, shape (n * m + k, n * m + k)
  189. Jacobian of the collocation system in a sparse form.
  190. References
  191. ----------
  192. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  193. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  194. Number 3, pp. 299-316, 2001.
  195. """
  196. df_dy = np.transpose(df_dy, (2, 0, 1))
  197. df_dy_middle = np.transpose(df_dy_middle, (2, 0, 1))
  198. h = h[:, np.newaxis, np.newaxis]
  199. dtype = df_dy.dtype
  200. # Computing diagonal n x n blocks.
  201. dPhi_dy_0 = np.empty((m - 1, n, n), dtype=dtype)
  202. dPhi_dy_0[:] = -np.identity(n)
  203. dPhi_dy_0 -= h / 6 * (df_dy[:-1] + 2 * df_dy_middle)
  204. T = stacked_matmul(df_dy_middle, df_dy[:-1])
  205. dPhi_dy_0 -= h**2 / 12 * T
  206. # Computing off-diagonal n x n blocks.
  207. dPhi_dy_1 = np.empty((m - 1, n, n), dtype=dtype)
  208. dPhi_dy_1[:] = np.identity(n)
  209. dPhi_dy_1 -= h / 6 * (df_dy[1:] + 2 * df_dy_middle)
  210. T = stacked_matmul(df_dy_middle, df_dy[1:])
  211. dPhi_dy_1 += h**2 / 12 * T
  212. values = np.hstack((dPhi_dy_0.ravel(), dPhi_dy_1.ravel(), dbc_dya.ravel(),
  213. dbc_dyb.ravel()))
  214. if k > 0:
  215. df_dp = np.transpose(df_dp, (2, 0, 1))
  216. df_dp_middle = np.transpose(df_dp_middle, (2, 0, 1))
  217. T = stacked_matmul(df_dy_middle, df_dp[:-1] - df_dp[1:])
  218. df_dp_middle += 0.125 * h * T
  219. dPhi_dp = -h/6 * (df_dp[:-1] + df_dp[1:] + 4 * df_dp_middle)
  220. values = np.hstack((values, dPhi_dp.ravel(), dbc_dp.ravel()))
  221. J = coo_matrix((values, (i_jac, j_jac)))
  222. return csc_matrix(J)
  223. def collocation_fun(fun, y, p, x, h):
  224. """Evaluate collocation residuals.
  225. This function lies in the core of the method. The solution is sought
  226. as a cubic C1 continuous spline with derivatives matching the ODE rhs
  227. at given nodes `x`. Collocation conditions are formed from the equality
  228. of the spline derivatives and rhs of the ODE system in the middle points
  229. between nodes.
  230. Such method is classified to Lobbato IIIA family in ODE literature.
  231. Refer to [1]_ for the formula and some discussion.
  232. Returns
  233. -------
  234. col_res : ndarray, shape (n, m - 1)
  235. Collocation residuals at the middle points of the mesh intervals.
  236. y_middle : ndarray, shape (n, m - 1)
  237. Values of the cubic spline evaluated at the middle points of the mesh
  238. intervals.
  239. f : ndarray, shape (n, m)
  240. RHS of the ODE system evaluated at the mesh nodes.
  241. f_middle : ndarray, shape (n, m - 1)
  242. RHS of the ODE system evaluated at the middle points of the mesh
  243. intervals (and using `y_middle`).
  244. References
  245. ----------
  246. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  247. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  248. Number 3, pp. 299-316, 2001.
  249. """
  250. f = fun(x, y, p)
  251. y_middle = (0.5 * (y[:, 1:] + y[:, :-1]) -
  252. 0.125 * h * (f[:, 1:] - f[:, :-1]))
  253. f_middle = fun(x[:-1] + 0.5 * h, y_middle, p)
  254. col_res = y[:, 1:] - y[:, :-1] - h / 6 * (f[:, :-1] + f[:, 1:] +
  255. 4 * f_middle)
  256. return col_res, y_middle, f, f_middle
  257. def prepare_sys(n, m, k, fun, bc, fun_jac, bc_jac, x, h):
  258. """Create the function and the Jacobian for the collocation system."""
  259. x_middle = x[:-1] + 0.5 * h
  260. i_jac, j_jac = compute_jac_indices(n, m, k)
  261. def col_fun(y, p):
  262. return collocation_fun(fun, y, p, x, h)
  263. def sys_jac(y, p, y_middle, f, f_middle, bc0):
  264. if fun_jac is None:
  265. df_dy, df_dp = estimate_fun_jac(fun, x, y, p, f)
  266. df_dy_middle, df_dp_middle = estimate_fun_jac(
  267. fun, x_middle, y_middle, p, f_middle)
  268. else:
  269. df_dy, df_dp = fun_jac(x, y, p)
  270. df_dy_middle, df_dp_middle = fun_jac(x_middle, y_middle, p)
  271. if bc_jac is None:
  272. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(bc, y[:, 0], y[:, -1],
  273. p, bc0)
  274. else:
  275. dbc_dya, dbc_dyb, dbc_dp = bc_jac(y[:, 0], y[:, -1], p)
  276. return construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy,
  277. df_dy_middle, df_dp, df_dp_middle, dbc_dya,
  278. dbc_dyb, dbc_dp)
  279. return col_fun, sys_jac
  280. def solve_newton(n, m, h, col_fun, bc, jac, y, p, B, bvp_tol, bc_tol):
  281. """Solve the nonlinear collocation system by a Newton method.
  282. This is a simple Newton method with a backtracking line search. As
  283. advised in [1]_, an affine-invariant criterion function F = ||J^-1 r||^2
  284. is used, where J is the Jacobian matrix at the current iteration and r is
  285. the vector or collocation residuals (values of the system lhs).
  286. The method alters between full Newton iterations and the fixed-Jacobian
  287. iterations based
  288. There are other tricks proposed in [1]_, but they are not used as they
  289. don't seem to improve anything significantly, and even break the
  290. convergence on some test problems I tried.
  291. All important parameters of the algorithm are defined inside the function.
  292. Parameters
  293. ----------
  294. n : int
  295. Number of equations in the ODE system.
  296. m : int
  297. Number of nodes in the mesh.
  298. h : ndarray, shape (m-1,)
  299. Mesh intervals.
  300. col_fun : callable
  301. Function computing collocation residuals.
  302. bc : callable
  303. Function computing boundary condition residuals.
  304. jac : callable
  305. Function computing the Jacobian of the whole system (including
  306. collocation and boundary condition residuals). It is supposed to
  307. return csc_matrix.
  308. y : ndarray, shape (n, m)
  309. Initial guess for the function values at the mesh nodes.
  310. p : ndarray, shape (k,)
  311. Initial guess for the unknown parameters.
  312. B : ndarray with shape (n, n) or None
  313. Matrix to force the S y(a) = 0 condition for a problems with the
  314. singular term. If None, the singular term is assumed to be absent.
  315. bvp_tol : float
  316. Tolerance to which we want to solve a BVP.
  317. bc_tol : float
  318. Tolerance to which we want to satisfy the boundary conditions.
  319. Returns
  320. -------
  321. y : ndarray, shape (n, m)
  322. Final iterate for the function values at the mesh nodes.
  323. p : ndarray, shape (k,)
  324. Final iterate for the unknown parameters.
  325. singular : bool
  326. True, if the LU decomposition failed because Jacobian turned out
  327. to be singular.
  328. References
  329. ----------
  330. .. [1] U. Ascher, R. Mattheij and R. Russell "Numerical Solution of
  331. Boundary Value Problems for Ordinary Differential Equations",
  332. Philidelphia, PA: Society for Industrial and Applied Mathematics,
  333. 1995.
  334. """
  335. # We know that the solution residuals at the middle points of the mesh
  336. # are connected with collocation residuals r_middle = 1.5 * col_res / h.
  337. # As our BVP solver tries to decrease relative residuals below a certain
  338. # tolerance, it seems reasonable to terminated Newton iterations by
  339. # comparison of r_middle / (1 + np.abs(f_middle)) with a certain threshold,
  340. # which we choose to be 1.5 orders lower than the BVP tolerance. We rewrite
  341. # the condition as col_res < tol_r * (1 + np.abs(f_middle)), then tol_r
  342. # should be computed as follows:
  343. tol_r = 2/3 * h * 5e-2 * bvp_tol
  344. # Maximum allowed number of Jacobian evaluation and factorization, in
  345. # other words, the maximum number of full Newton iterations. A small value
  346. # is recommended in the literature.
  347. max_njev = 4
  348. # Maximum number of iterations, considering that some of them can be
  349. # performed with the fixed Jacobian. In theory, such iterations are cheap,
  350. # but it's not that simple in Python.
  351. max_iter = 8
  352. # Minimum relative improvement of the criterion function to accept the
  353. # step (Armijo constant).
  354. sigma = 0.2
  355. # Step size decrease factor for backtracking.
  356. tau = 0.5
  357. # Maximum number of backtracking steps, the minimum step is then
  358. # tau ** n_trial.
  359. n_trial = 4
  360. col_res, y_middle, f, f_middle = col_fun(y, p)
  361. bc_res = bc(y[:, 0], y[:, -1], p)
  362. res = np.hstack((col_res.ravel(order='F'), bc_res))
  363. njev = 0
  364. singular = False
  365. recompute_jac = True
  366. for iteration in range(max_iter):
  367. if recompute_jac:
  368. J = jac(y, p, y_middle, f, f_middle, bc_res)
  369. njev += 1
  370. try:
  371. LU = splu(J)
  372. except RuntimeError:
  373. singular = True
  374. break
  375. step = LU.solve(res)
  376. cost = np.dot(step, step)
  377. y_step = step[:m * n].reshape((n, m), order='F')
  378. p_step = step[m * n:]
  379. alpha = 1
  380. for trial in range(n_trial + 1):
  381. y_new = y - alpha * y_step
  382. if B is not None:
  383. y_new[:, 0] = np.dot(B, y_new[:, 0])
  384. p_new = p - alpha * p_step
  385. col_res, y_middle, f, f_middle = col_fun(y_new, p_new)
  386. bc_res = bc(y_new[:, 0], y_new[:, -1], p_new)
  387. res = np.hstack((col_res.ravel(order='F'), bc_res))
  388. step_new = LU.solve(res)
  389. cost_new = np.dot(step_new, step_new)
  390. if cost_new < (1 - 2 * alpha * sigma) * cost:
  391. break
  392. if trial < n_trial:
  393. alpha *= tau
  394. y = y_new
  395. p = p_new
  396. if njev == max_njev:
  397. break
  398. if (np.all(np.abs(col_res) < tol_r * (1 + np.abs(f_middle))) and
  399. np.all(np.abs(bc_res) < bc_tol)):
  400. break
  401. # If the full step was taken, then we are going to continue with
  402. # the same Jacobian. This is the approach of BVP_SOLVER.
  403. if alpha == 1:
  404. step = step_new
  405. cost = cost_new
  406. recompute_jac = False
  407. else:
  408. recompute_jac = True
  409. return y, p, singular
  410. def print_iteration_header():
  411. print(f"{'Iteration':^15}{'Max residual':^15}{'Max BC residual':^15}"
  412. f"{'Total nodes':^15}{'Nodes added':^15}")
  413. def print_iteration_progress(iteration, residual, bc_residual, total_nodes,
  414. nodes_added):
  415. print(f"{iteration:^15}{residual:^15.2e}{bc_residual:^15.2e}"
  416. f"{total_nodes:^15}{nodes_added:^15}")
  417. class BVPResult(OptimizeResult):
  418. pass
  419. TERMINATION_MESSAGES = {
  420. 0: "The algorithm converged to the desired accuracy.",
  421. 1: "The maximum number of mesh nodes is exceeded.",
  422. 2: "A singular Jacobian encountered when solving the collocation system.",
  423. 3: "The solver was unable to satisfy boundary conditions tolerance on iteration 10."
  424. }
  425. def estimate_rms_residuals(fun, sol, x, h, p, r_middle, f_middle):
  426. """Estimate rms values of collocation residuals using Lobatto quadrature.
  427. The residuals are defined as the difference between the derivatives of
  428. our solution and rhs of the ODE system. We use relative residuals, i.e.,
  429. normalized by 1 + np.abs(f). RMS values are computed as sqrt from the
  430. normalized integrals of the squared relative residuals over each interval.
  431. Integrals are estimated using 5-point Lobatto quadrature [1]_, we use the
  432. fact that residuals at the mesh nodes are identically zero.
  433. In [2] they don't normalize integrals by interval lengths, which gives
  434. a higher rate of convergence of the residuals by the factor of h**0.5.
  435. I chose to do such normalization for an ease of interpretation of return
  436. values as RMS estimates.
  437. Returns
  438. -------
  439. rms_res : ndarray, shape (m - 1,)
  440. Estimated rms values of the relative residuals over each interval.
  441. References
  442. ----------
  443. .. [1] http://mathworld.wolfram.com/LobattoQuadrature.html
  444. .. [2] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  445. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  446. Number 3, pp. 299-316, 2001.
  447. """
  448. x_middle = x[:-1] + 0.5 * h
  449. s = 0.5 * h * (3/7)**0.5
  450. x1 = x_middle + s
  451. x2 = x_middle - s
  452. y1 = sol(x1)
  453. y2 = sol(x2)
  454. y1_prime = sol(x1, 1)
  455. y2_prime = sol(x2, 1)
  456. f1 = fun(x1, y1, p)
  457. f2 = fun(x2, y2, p)
  458. r1 = y1_prime - f1
  459. r2 = y2_prime - f2
  460. r_middle /= 1 + np.abs(f_middle)
  461. r1 /= 1 + np.abs(f1)
  462. r2 /= 1 + np.abs(f2)
  463. r1 = np.sum(np.real(r1 * np.conj(r1)), axis=0)
  464. r2 = np.sum(np.real(r2 * np.conj(r2)), axis=0)
  465. r_middle = np.sum(np.real(r_middle * np.conj(r_middle)), axis=0)
  466. return (0.5 * (32 / 45 * r_middle + 49 / 90 * (r1 + r2))) ** 0.5
  467. def create_spline(y, yp, x, h):
  468. """Create a cubic spline given values and derivatives.
  469. Formulas for the coefficients are taken from interpolate.CubicSpline.
  470. Returns
  471. -------
  472. sol : PPoly
  473. Constructed spline as a PPoly instance.
  474. """
  475. from scipy.interpolate import PPoly
  476. n, m = y.shape
  477. c = np.empty((4, n, m - 1), dtype=y.dtype)
  478. slope = (y[:, 1:] - y[:, :-1]) / h
  479. t = (yp[:, :-1] + yp[:, 1:] - 2 * slope) / h
  480. c[0] = t / h
  481. c[1] = (slope - yp[:, :-1]) / h - t
  482. c[2] = yp[:, :-1]
  483. c[3] = y[:, :-1]
  484. c = np.moveaxis(c, 1, 0)
  485. return PPoly(c, x, extrapolate=True, axis=1)
  486. def modify_mesh(x, insert_1, insert_2):
  487. """Insert nodes into a mesh.
  488. Nodes removal logic is not established, its impact on the solver is
  489. presumably negligible. So, only insertion is done in this function.
  490. Parameters
  491. ----------
  492. x : ndarray, shape (m,)
  493. Mesh nodes.
  494. insert_1 : ndarray
  495. Intervals to each insert 1 new node in the middle.
  496. insert_2 : ndarray
  497. Intervals to each insert 2 new nodes, such that divide an interval
  498. into 3 equal parts.
  499. Returns
  500. -------
  501. x_new : ndarray
  502. New mesh nodes.
  503. Notes
  504. -----
  505. `insert_1` and `insert_2` should not have common values.
  506. """
  507. # Because np.insert implementation apparently varies with a version of
  508. # NumPy, we use a simple and reliable approach with sorting.
  509. return np.sort(np.hstack((
  510. x,
  511. 0.5 * (x[insert_1] + x[insert_1 + 1]),
  512. (2 * x[insert_2] + x[insert_2 + 1]) / 3,
  513. (x[insert_2] + 2 * x[insert_2 + 1]) / 3
  514. )))
  515. def wrap_functions(fun, bc, fun_jac, bc_jac, k, a, S, D, dtype):
  516. """Wrap functions for unified usage in the solver."""
  517. if fun_jac is None:
  518. fun_jac_wrapped = None
  519. if bc_jac is None:
  520. bc_jac_wrapped = None
  521. if k == 0:
  522. def fun_p(x, y, _):
  523. return np.asarray(fun(x, y), dtype)
  524. def bc_wrapped(ya, yb, _):
  525. return np.asarray(bc(ya, yb), dtype)
  526. if fun_jac is not None:
  527. def fun_jac_p(x, y, _):
  528. return np.asarray(fun_jac(x, y), dtype), None
  529. if bc_jac is not None:
  530. def bc_jac_wrapped(ya, yb, _):
  531. dbc_dya, dbc_dyb = bc_jac(ya, yb)
  532. return (np.asarray(dbc_dya, dtype),
  533. np.asarray(dbc_dyb, dtype), None)
  534. else:
  535. def fun_p(x, y, p):
  536. return np.asarray(fun(x, y, p), dtype)
  537. def bc_wrapped(x, y, p):
  538. return np.asarray(bc(x, y, p), dtype)
  539. if fun_jac is not None:
  540. def fun_jac_p(x, y, p):
  541. df_dy, df_dp = fun_jac(x, y, p)
  542. return np.asarray(df_dy, dtype), np.asarray(df_dp, dtype)
  543. if bc_jac is not None:
  544. def bc_jac_wrapped(ya, yb, p):
  545. dbc_dya, dbc_dyb, dbc_dp = bc_jac(ya, yb, p)
  546. return (np.asarray(dbc_dya, dtype), np.asarray(dbc_dyb, dtype),
  547. np.asarray(dbc_dp, dtype))
  548. if S is None:
  549. fun_wrapped = fun_p
  550. else:
  551. def fun_wrapped(x, y, p):
  552. f = fun_p(x, y, p)
  553. if x[0] == a:
  554. f[:, 0] = np.dot(D, f[:, 0])
  555. f[:, 1:] += np.dot(S, y[:, 1:]) / (x[1:] - a)
  556. else:
  557. f += np.dot(S, y) / (x - a)
  558. return f
  559. if fun_jac is not None:
  560. if S is None:
  561. fun_jac_wrapped = fun_jac_p
  562. else:
  563. Sr = S[:, :, np.newaxis]
  564. def fun_jac_wrapped(x, y, p):
  565. df_dy, df_dp = fun_jac_p(x, y, p)
  566. if x[0] == a:
  567. df_dy[:, :, 0] = np.dot(D, df_dy[:, :, 0])
  568. df_dy[:, :, 1:] += Sr / (x[1:] - a)
  569. else:
  570. df_dy += Sr / (x - a)
  571. return df_dy, df_dp
  572. return fun_wrapped, bc_wrapped, fun_jac_wrapped, bc_jac_wrapped
  573. @xp_capabilities(np_only=True)
  574. def solve_bvp(fun, bc, x, y, p=None, S=None, fun_jac=None, bc_jac=None,
  575. tol=1e-3, max_nodes=1000, verbose=0, bc_tol=None):
  576. """Solve a boundary value problem for a system of ODEs.
  577. This function numerically solves a first order system of ODEs subject to
  578. two-point boundary conditions::
  579. dy / dx = f(x, y, p) + S * y / (x - a), a <= x <= b
  580. bc(y(a), y(b), p) = 0
  581. Here x is a 1-D independent variable, y(x) is an n-D
  582. vector-valued function and p is a k-D vector of unknown
  583. parameters which is to be found along with y(x). For the problem to be
  584. determined, there must be n + k boundary conditions, i.e., bc must be an
  585. (n + k)-D function.
  586. The last singular term on the right-hand side of the system is optional.
  587. It is defined by an n-by-n matrix S, such that the solution must satisfy
  588. S y(a) = 0. This condition will be forced during iterations, so it must not
  589. contradict boundary conditions. See [2]_ for the explanation how this term
  590. is handled when solving BVPs numerically.
  591. Problems in a complex domain can be solved as well. In this case, y and p
  592. are considered to be complex, and f and bc are assumed to be complex-valued
  593. functions, but x stays real. Note that f and bc must be complex
  594. differentiable (satisfy Cauchy-Riemann equations [4]_), otherwise you
  595. should rewrite your problem for real and imaginary parts separately. To
  596. solve a problem in a complex domain, pass an initial guess for y with a
  597. complex data type (see below).
  598. Parameters
  599. ----------
  600. fun : callable
  601. Right-hand side of the system. The calling signature is ``fun(x, y)``,
  602. or ``fun(x, y, p)`` if parameters are present. All arguments are
  603. ndarray: ``x`` with shape (m,), ``y`` with shape (n, m), meaning that
  604. ``y[:, i]`` corresponds to ``x[i]``, and ``p`` with shape (k,). The
  605. return value must be an array with shape (n, m) and with the same
  606. layout as ``y``.
  607. bc : callable
  608. Function evaluating residuals of the boundary conditions. The calling
  609. signature is ``bc(ya, yb)``, or ``bc(ya, yb, p)`` if parameters are
  610. present. All arguments are ndarray: ``ya`` and ``yb`` with shape (n,),
  611. and ``p`` with shape (k,). The return value must be an array with
  612. shape (n + k,).
  613. x : array_like, shape (m,)
  614. Initial mesh. Must be a strictly increasing sequence of real numbers
  615. with ``x[0]=a`` and ``x[-1]=b``.
  616. y : array_like, shape (n, m)
  617. Initial guess for the function values at the mesh nodes, ith column
  618. corresponds to ``x[i]``. For problems in a complex domain pass `y`
  619. with a complex data type (even if the initial guess is purely real).
  620. p : array_like with shape (k,) or None, optional
  621. Initial guess for the unknown parameters. If None (default), it is
  622. assumed that the problem doesn't depend on any parameters.
  623. S : array_like with shape (n, n) or None
  624. Matrix defining the singular term. If None (default), the problem is
  625. solved without the singular term.
  626. fun_jac : callable or None, optional
  627. Function computing derivatives of f with respect to y and p. The
  628. calling signature is ``fun_jac(x, y)``, or ``fun_jac(x, y, p)`` if
  629. parameters are present. The return must contain 1 or 2 elements in the
  630. following order:
  631. * df_dy : array_like with shape (n, n, m), where an element
  632. (i, j, q) equals to d f_i(x_q, y_q, p) / d (y_q)_j.
  633. * df_dp : array_like with shape (n, k, m), where an element
  634. (i, j, q) equals to d f_i(x_q, y_q, p) / d p_j.
  635. Here q numbers nodes at which x and y are defined, whereas i and j
  636. number vector components. If the problem is solved without unknown
  637. parameters, df_dp should not be returned.
  638. If `fun_jac` is None (default), the derivatives will be estimated
  639. by the forward finite differences.
  640. bc_jac : callable or None, optional
  641. Function computing derivatives of bc with respect to ya, yb, and p.
  642. The calling signature is ``bc_jac(ya, yb)``, or ``bc_jac(ya, yb, p)``
  643. if parameters are present. The return must contain 2 or 3 elements in
  644. the following order:
  645. * dbc_dya : array_like with shape (n, n), where an element (i, j)
  646. equals to d bc_i(ya, yb, p) / d ya_j.
  647. * dbc_dyb : array_like with shape (n, n), where an element (i, j)
  648. equals to d bc_i(ya, yb, p) / d yb_j.
  649. * dbc_dp : array_like with shape (n, k), where an element (i, j)
  650. equals to d bc_i(ya, yb, p) / d p_j.
  651. If the problem is solved without unknown parameters, dbc_dp should not
  652. be returned.
  653. If `bc_jac` is None (default), the derivatives will be estimated by
  654. the forward finite differences.
  655. tol : float, optional
  656. Desired tolerance of the solution. If we define ``r = y' - f(x, y)``,
  657. where y is the found solution, then the solver tries to achieve on each
  658. mesh interval ``norm(r / (1 + abs(f)) < tol``, where ``norm`` is
  659. estimated in a root mean squared sense (using a numerical quadrature
  660. formula). Default is 1e-3.
  661. max_nodes : int, optional
  662. Maximum allowed number of the mesh nodes. If exceeded, the algorithm
  663. terminates. Default is 1000.
  664. verbose : {0, 1, 2}, optional
  665. Level of algorithm's verbosity:
  666. * 0 (default) : work silently.
  667. * 1 : display a termination report.
  668. * 2 : display progress during iterations.
  669. bc_tol : float, optional
  670. Desired absolute tolerance for the boundary condition residuals: `bc`
  671. value should satisfy ``abs(bc) < bc_tol`` component-wise.
  672. Equals to `tol` by default. Up to 10 iterations are allowed to achieve this
  673. tolerance.
  674. Returns
  675. -------
  676. Bunch object with the following fields defined:
  677. sol : PPoly
  678. Found solution for y as `scipy.interpolate.PPoly` instance, a C1
  679. continuous cubic spline.
  680. p : ndarray or None, shape (k,)
  681. Found parameters. None, if the parameters were not present in the
  682. problem.
  683. x : ndarray, shape (m,)
  684. Nodes of the final mesh.
  685. y : ndarray, shape (n, m)
  686. Solution values at the mesh nodes.
  687. yp : ndarray, shape (n, m)
  688. Solution derivatives at the mesh nodes.
  689. rms_residuals : ndarray, shape (m - 1,)
  690. RMS values of the relative residuals over each mesh interval (see the
  691. description of `tol` parameter).
  692. niter : int
  693. Number of completed iterations.
  694. status : int
  695. Reason for algorithm termination:
  696. * 0: The algorithm converged to the desired accuracy.
  697. * 1: The maximum number of mesh nodes is exceeded.
  698. * 2: A singular Jacobian encountered when solving the collocation
  699. system.
  700. message : string
  701. Verbal description of the termination reason.
  702. success : bool
  703. True if the algorithm converged to the desired accuracy (``status=0``).
  704. Notes
  705. -----
  706. This function implements a 4th order collocation algorithm with the
  707. control of residuals similar to [1]_. A collocation system is solved
  708. by a damped Newton method with an affine-invariant criterion function as
  709. described in [3]_.
  710. Note that in [1]_ integral residuals are defined without normalization
  711. by interval lengths. So, their definition is different by a multiplier of
  712. h**0.5 (h is an interval length) from the definition used here.
  713. .. versionadded:: 0.18.0
  714. References
  715. ----------
  716. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  717. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  718. Number 3, pp. 299-316, 2001.
  719. .. [2] L.F. Shampine, P. H. Muir and H. Xu, "A User-Friendly Fortran BVP
  720. Solver", J. Numer. Anal., Ind. Appl. Math. (JNAIAM), Vol. 1,
  721. Number 2, pp. 201-217, 2006.
  722. .. [3] U. Ascher, R. Mattheij and R. Russell "Numerical Solution of
  723. Boundary Value Problems for Ordinary Differential Equations",
  724. Philidelphia, PA: Society for Industrial and Applied Mathematics,
  725. 1995.
  726. :doi:`10.1137/1.9781611971231`
  727. .. [4] `Cauchy-Riemann equations
  728. <https://en.wikipedia.org/wiki/Cauchy-Riemann_equations>`_ on
  729. Wikipedia.
  730. Examples
  731. --------
  732. In the first example, we solve Bratu's problem::
  733. y'' + k * exp(y) = 0
  734. y(0) = y(1) = 0
  735. for k = 1.
  736. We rewrite the equation as a first-order system and implement its
  737. right-hand side evaluation::
  738. y1' = y2
  739. y2' = -exp(y1)
  740. >>> import numpy as np
  741. >>> def fun(x, y):
  742. ... return np.vstack((y[1], -np.exp(y[0])))
  743. Implement evaluation of the boundary condition residuals:
  744. >>> def bc(ya, yb):
  745. ... return np.array([ya[0], yb[0]])
  746. Define the initial mesh with 5 nodes:
  747. >>> x = np.linspace(0, 1, 5)
  748. This problem is known to have two solutions. To obtain both of them, we
  749. use two different initial guesses for y. We denote them by subscripts
  750. a and b.
  751. >>> y_a = np.zeros((2, x.size))
  752. >>> y_b = np.zeros((2, x.size))
  753. >>> y_b[0] = 3
  754. Now we are ready to run the solver.
  755. >>> from scipy.integrate import solve_bvp
  756. >>> res_a = solve_bvp(fun, bc, x, y_a)
  757. >>> res_b = solve_bvp(fun, bc, x, y_b)
  758. Let's plot the two found solutions. We take an advantage of having the
  759. solution in a spline form to produce a smooth plot.
  760. >>> x_plot = np.linspace(0, 1, 100)
  761. >>> y_plot_a = res_a.sol(x_plot)[0]
  762. >>> y_plot_b = res_b.sol(x_plot)[0]
  763. >>> import matplotlib.pyplot as plt
  764. >>> plt.plot(x_plot, y_plot_a, label='y_a')
  765. >>> plt.plot(x_plot, y_plot_b, label='y_b')
  766. >>> plt.legend()
  767. >>> plt.xlabel("x")
  768. >>> plt.ylabel("y")
  769. >>> plt.show()
  770. We see that the two solutions have similar shape, but differ in scale
  771. significantly.
  772. In the second example, we solve a simple Sturm-Liouville problem::
  773. y'' + k**2 * y = 0
  774. y(0) = y(1) = 0
  775. It is known that a non-trivial solution y = A * sin(k * x) is possible for
  776. k = pi * n, where n is an integer. To establish the normalization constant
  777. A = 1 we add a boundary condition::
  778. y'(0) = k
  779. Again, we rewrite our equation as a first-order system and implement its
  780. right-hand side evaluation::
  781. y1' = y2
  782. y2' = -k**2 * y1
  783. >>> def fun(x, y, p):
  784. ... k = p[0]
  785. ... return np.vstack((y[1], -k**2 * y[0]))
  786. Note that parameters p are passed as a vector (with one element in our
  787. case).
  788. Implement the boundary conditions:
  789. >>> def bc(ya, yb, p):
  790. ... k = p[0]
  791. ... return np.array([ya[0], yb[0], ya[1] - k])
  792. Set up the initial mesh and guess for y. We aim to find the solution for
  793. k = 2 * pi, to achieve that we set values of y to approximately follow
  794. sin(2 * pi * x):
  795. >>> x = np.linspace(0, 1, 5)
  796. >>> y = np.zeros((2, x.size))
  797. >>> y[0, 1] = 1
  798. >>> y[0, 3] = -1
  799. Run the solver with 6 as an initial guess for k.
  800. >>> sol = solve_bvp(fun, bc, x, y, p=[6])
  801. We see that the found k is approximately correct:
  802. >>> sol.p[0]
  803. 6.28329460046
  804. And, finally, plot the solution to see the anticipated sinusoid:
  805. >>> x_plot = np.linspace(0, 1, 100)
  806. >>> y_plot = sol.sol(x_plot)[0]
  807. >>> plt.plot(x_plot, y_plot)
  808. >>> plt.xlabel("x")
  809. >>> plt.ylabel("y")
  810. >>> plt.show()
  811. """
  812. x = np.asarray(x, dtype=float)
  813. if x.ndim != 1:
  814. raise ValueError("`x` must be 1 dimensional.")
  815. h = np.diff(x)
  816. if np.any(h <= 0):
  817. raise ValueError("`x` must be strictly increasing.")
  818. a = x[0]
  819. y = np.asarray(y)
  820. if np.issubdtype(y.dtype, np.complexfloating):
  821. dtype = complex
  822. else:
  823. dtype = float
  824. y = y.astype(dtype, copy=False)
  825. if y.ndim != 2:
  826. raise ValueError("`y` must be 2 dimensional.")
  827. if y.shape[1] != x.shape[0]:
  828. raise ValueError(f"`y` is expected to have {x.shape[0]} columns, but actually "
  829. f"has {y.shape[1]}.")
  830. if p is None:
  831. p = np.array([])
  832. else:
  833. p = np.asarray(p, dtype=dtype)
  834. if p.ndim != 1:
  835. raise ValueError("`p` must be 1 dimensional.")
  836. if tol < 100 * EPS:
  837. warn(f"`tol` is too low, setting to {100 * EPS:.2e}", stacklevel=2)
  838. tol = 100 * EPS
  839. if verbose not in [0, 1, 2]:
  840. raise ValueError("`verbose` must be in [0, 1, 2].")
  841. n = y.shape[0]
  842. k = p.shape[0]
  843. if S is not None:
  844. S = np.asarray(S, dtype=dtype)
  845. if S.shape != (n, n):
  846. raise ValueError(f"`S` is expected to have shape {(n, n)}, "
  847. f"but actually has {S.shape}")
  848. # Compute I - S^+ S to impose necessary boundary conditions.
  849. B = np.identity(n) - np.dot(pinv(S), S)
  850. y[:, 0] = np.dot(B, y[:, 0])
  851. # Compute (I - S)^+ to correct derivatives at x=a.
  852. D = pinv(np.identity(n) - S)
  853. else:
  854. B = None
  855. D = None
  856. if bc_tol is None:
  857. bc_tol = tol
  858. # Maximum number of iterations
  859. max_iteration = 10
  860. fun_wrapped, bc_wrapped, fun_jac_wrapped, bc_jac_wrapped = wrap_functions(
  861. fun, bc, fun_jac, bc_jac, k, a, S, D, dtype)
  862. f = fun_wrapped(x, y, p)
  863. if f.shape != y.shape:
  864. raise ValueError(f"`fun` return is expected to have shape {y.shape}, "
  865. f"but actually has {f.shape}.")
  866. bc_res = bc_wrapped(y[:, 0], y[:, -1], p)
  867. if bc_res.shape != (n + k,):
  868. raise ValueError(f"`bc` return is expected to have shape {(n + k,)}, "
  869. f"but actually has {bc_res.shape}.")
  870. status = 0
  871. iteration = 0
  872. if verbose == 2:
  873. print_iteration_header()
  874. while True:
  875. m = x.shape[0]
  876. col_fun, jac_sys = prepare_sys(n, m, k, fun_wrapped, bc_wrapped,
  877. fun_jac_wrapped, bc_jac_wrapped, x, h)
  878. y, p, singular = solve_newton(n, m, h, col_fun, bc_wrapped, jac_sys,
  879. y, p, B, tol, bc_tol)
  880. iteration += 1
  881. col_res, y_middle, f, f_middle = collocation_fun(fun_wrapped, y,
  882. p, x, h)
  883. bc_res = bc_wrapped(y[:, 0], y[:, -1], p)
  884. max_bc_res = np.max(abs(bc_res))
  885. # This relation is not trivial, but can be verified.
  886. r_middle = 1.5 * col_res / h
  887. sol = create_spline(y, f, x, h)
  888. rms_res = estimate_rms_residuals(fun_wrapped, sol, x, h, p,
  889. r_middle, f_middle)
  890. max_rms_res = np.max(rms_res)
  891. if singular:
  892. status = 2
  893. break
  894. insert_1, = np.nonzero((rms_res > tol) & (rms_res < 100 * tol))
  895. insert_2, = np.nonzero(rms_res >= 100 * tol)
  896. nodes_added = insert_1.shape[0] + 2 * insert_2.shape[0]
  897. if m + nodes_added > max_nodes:
  898. status = 1
  899. if verbose == 2:
  900. nodes_added = f"({nodes_added})"
  901. print_iteration_progress(iteration, max_rms_res, max_bc_res,
  902. m, nodes_added)
  903. break
  904. if verbose == 2:
  905. print_iteration_progress(iteration, max_rms_res, max_bc_res, m,
  906. nodes_added)
  907. if nodes_added > 0:
  908. x = modify_mesh(x, insert_1, insert_2)
  909. h = np.diff(x)
  910. y = sol(x)
  911. elif max_bc_res <= bc_tol:
  912. status = 0
  913. break
  914. elif iteration >= max_iteration:
  915. status = 3
  916. break
  917. if verbose > 0:
  918. if status == 0:
  919. print(f"Solved in {iteration} iterations, number of nodes {x.shape[0]}. \n"
  920. f"Maximum relative residual: {max_rms_res:.2e} \n"
  921. f"Maximum boundary residual: {max_bc_res:.2e}")
  922. elif status == 1:
  923. print(f"Number of nodes is exceeded after iteration {iteration}. \n"
  924. f"Maximum relative residual: {max_rms_res:.2e} \n"
  925. f"Maximum boundary residual: {max_bc_res:.2e}")
  926. elif status == 2:
  927. print("Singular Jacobian encountered when solving the collocation "
  928. f"system on iteration {iteration}. \n"
  929. f"Maximum relative residual: {max_rms_res:.2e} \n"
  930. f"Maximum boundary residual: {max_bc_res:.2e}")
  931. elif status == 3:
  932. print("The solver was unable to satisfy boundary conditions "
  933. f"tolerance on iteration {iteration}. \n"
  934. f"Maximum relative residual: {max_rms_res:.2e} \n"
  935. f"Maximum boundary residual: {max_bc_res:.2e}")
  936. if p.size == 0:
  937. p = None
  938. return BVPResult(sol=sol, p=p, x=x, y=y, yp=f, rms_residuals=rms_res,
  939. niter=iteration, status=status,
  940. message=TERMINATION_MESSAGES[status], success=status == 0)