_bvp.py 40 KB

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