_bvp.py 39 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134
  1. """Boundary value problem solver."""
  2. from __future__ import division, print_function, absolute_import
  3. from warnings import warn
  4. import numpy as np
  5. from numpy.linalg import norm, pinv
  6. from scipy.sparse import coo_matrix, csc_matrix
  7. from scipy.sparse.linalg import splu
  8. from scipy.optimize import OptimizeResult
  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. In our case a[i, :, :] and b[i, :, :] are always square.
  117. """
  118. # Empirical optimization. Use outer Python loop and BLAS for large
  119. # matrices, otherwise use a single einsum call.
  120. if a.shape[1] > 50:
  121. out = np.empty_like(a)
  122. for i in range(a.shape[0]):
  123. out[i] = np.dot(a[i], b[i])
  124. return out
  125. else:
  126. return np.einsum('...ij,...jk->...ik', a, b)
  127. def construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy, df_dy_middle, df_dp,
  128. df_dp_middle, dbc_dya, dbc_dyb, dbc_dp):
  129. """Construct the Jacobian of the collocation system.
  130. There are n * m + k functions: m - 1 collocations residuals, each
  131. containing n components, followed by n + k boundary condition residuals.
  132. There are n * m + k variables: m vectors of y, each containing n
  133. components, followed by k values of vector p.
  134. For example, let m = 4, n = 2 and k = 1, then the Jacobian will have
  135. the following sparsity structure:
  136. 1 1 2 2 0 0 0 0 5
  137. 1 1 2 2 0 0 0 0 5
  138. 0 0 1 1 2 2 0 0 5
  139. 0 0 1 1 2 2 0 0 5
  140. 0 0 0 0 1 1 2 2 5
  141. 0 0 0 0 1 1 2 2 5
  142. 3 3 0 0 0 0 4 4 6
  143. 3 3 0 0 0 0 4 4 6
  144. 3 3 0 0 0 0 4 4 6
  145. Zeros denote identically zero values, other values denote different kinds
  146. of blocks in the matrix (see below). The blank row indicates the separation
  147. of collocation residuals from boundary conditions. And the blank column
  148. indicates the separation of y values from p values.
  149. Refer to [1]_ (p. 306) for the formula of n x n blocks for derivatives
  150. of collocation residuals with respect to y.
  151. Parameters
  152. ----------
  153. n : int
  154. Number of equations in the ODE system.
  155. m : int
  156. Number of nodes in the mesh.
  157. k : int
  158. Number of the unknown parameters.
  159. i_jac, j_jac : ndarray
  160. Row and column indices returned by `compute_jac_indices`. They
  161. represent different blocks in the Jacobian matrix in the following
  162. order (see the scheme above):
  163. * 1: m - 1 diagonal n x n blocks for the collocation residuals.
  164. * 2: m - 1 off-diagonal n x n blocks for the collocation residuals.
  165. * 3 : (n + k) x n block for the dependency of the boundary
  166. conditions on ya.
  167. * 4: (n + k) x n block for the dependency of the boundary
  168. conditions on yb.
  169. * 5: (m - 1) * n x k block for the dependency of the collocation
  170. residuals on p.
  171. * 6: (n + k) x k block for the dependency of the boundary
  172. conditions on p.
  173. df_dy : ndarray, shape (n, n, m)
  174. Jacobian of f with respect to y computed at the mesh nodes.
  175. df_dy_middle : ndarray, shape (n, n, m - 1)
  176. Jacobian of f with respect to y computed at the middle between the
  177. mesh nodes.
  178. df_dp : ndarray with shape (n, k, m) or None
  179. Jacobian of f with respect to p computed at the mesh nodes.
  180. df_dp_middle: ndarray with shape (n, k, m - 1) or None
  181. Jacobian of f with respect to p computed at the middle between the
  182. mesh nodes.
  183. dbc_dya, dbc_dyb : ndarray, shape (n, n)
  184. Jacobian of bc with respect to ya and yb.
  185. dbc_dp: ndarray with shape (n, k) or None
  186. Jacobian of bc with respect to p.
  187. Returns
  188. -------
  189. J : csc_matrix, shape (n * m + k, n * m + k)
  190. Jacobian of the collocation system in a sparse form.
  191. References
  192. ----------
  193. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  194. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  195. Number 3, pp. 299-316, 2001.
  196. """
  197. df_dy = np.transpose(df_dy, (2, 0, 1))
  198. df_dy_middle = np.transpose(df_dy_middle, (2, 0, 1))
  199. h = h[:, np.newaxis, np.newaxis]
  200. dtype = df_dy.dtype
  201. # Computing diagonal n x n blocks.
  202. dPhi_dy_0 = np.empty((m - 1, n, n), dtype=dtype)
  203. dPhi_dy_0[:] = -np.identity(n)
  204. dPhi_dy_0 -= h / 6 * (df_dy[:-1] + 2 * df_dy_middle)
  205. T = stacked_matmul(df_dy_middle, df_dy[:-1])
  206. dPhi_dy_0 -= h**2 / 12 * T
  207. # Computing off-diagonal n x n blocks.
  208. dPhi_dy_1 = np.empty((m - 1, n, n), dtype=dtype)
  209. dPhi_dy_1[:] = np.identity(n)
  210. dPhi_dy_1 -= h / 6 * (df_dy[1:] + 2 * df_dy_middle)
  211. T = stacked_matmul(df_dy_middle, df_dy[1:])
  212. dPhi_dy_1 += h**2 / 12 * T
  213. values = np.hstack((dPhi_dy_0.ravel(), dPhi_dy_1.ravel(), dbc_dya.ravel(),
  214. dbc_dyb.ravel()))
  215. if k > 0:
  216. df_dp = np.transpose(df_dp, (2, 0, 1))
  217. df_dp_middle = np.transpose(df_dp_middle, (2, 0, 1))
  218. T = stacked_matmul(df_dy_middle, df_dp[:-1] - df_dp[1:])
  219. df_dp_middle += 0.125 * h * T
  220. dPhi_dp = -h/6 * (df_dp[:-1] + df_dp[1:] + 4 * df_dp_middle)
  221. values = np.hstack((values, dPhi_dp.ravel(), dbc_dp.ravel()))
  222. J = coo_matrix((values, (i_jac, j_jac)))
  223. return csc_matrix(J)
  224. def collocation_fun(fun, y, p, x, h):
  225. """Evaluate collocation residuals.
  226. This function lies in the core of the method. The solution is sought
  227. as a cubic C1 continuous spline with derivatives matching the ODE rhs
  228. at given nodes `x`. Collocation conditions are formed from the equality
  229. of the spline derivatives and rhs of the ODE system in the middle points
  230. between nodes.
  231. Such method is classified to Lobbato IIIA family in ODE literature.
  232. Refer to [1]_ for the formula and some discussion.
  233. Returns
  234. -------
  235. col_res : ndarray, shape (n, m - 1)
  236. Collocation residuals at the middle points of the mesh intervals.
  237. y_middle : ndarray, shape (n, m - 1)
  238. Values of the cubic spline evaluated at the middle points of the mesh
  239. intervals.
  240. f : ndarray, shape (n, m)
  241. RHS of the ODE system evaluated at the mesh nodes.
  242. f_middle : ndarray, shape (n, m - 1)
  243. RHS of the ODE system evaluated at the middle points of the mesh
  244. intervals (and using `y_middle`).
  245. References
  246. ----------
  247. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  248. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  249. Number 3, pp. 299-316, 2001.
  250. """
  251. f = fun(x, y, p)
  252. y_middle = (0.5 * (y[:, 1:] + y[:, :-1]) -
  253. 0.125 * h * (f[:, 1:] - f[:, :-1]))
  254. f_middle = fun(x[:-1] + 0.5 * h, y_middle, p)
  255. col_res = y[:, 1:] - y[:, :-1] - h / 6 * (f[:, :-1] + f[:, 1:] +
  256. 4 * f_middle)
  257. return col_res, y_middle, f, f_middle
  258. def prepare_sys(n, m, k, fun, bc, fun_jac, bc_jac, x, h):
  259. """Create the function and the Jacobian for the collocation system."""
  260. x_middle = x[:-1] + 0.5 * h
  261. i_jac, j_jac = compute_jac_indices(n, m, k)
  262. def col_fun(y, p):
  263. return collocation_fun(fun, y, p, x, h)
  264. def sys_jac(y, p, y_middle, f, f_middle, bc0):
  265. if fun_jac is None:
  266. df_dy, df_dp = estimate_fun_jac(fun, x, y, p, f)
  267. df_dy_middle, df_dp_middle = estimate_fun_jac(
  268. fun, x_middle, y_middle, p, f_middle)
  269. else:
  270. df_dy, df_dp = fun_jac(x, y, p)
  271. df_dy_middle, df_dp_middle = fun_jac(x_middle, y_middle, p)
  272. if bc_jac is None:
  273. dbc_dya, dbc_dyb, dbc_dp = estimate_bc_jac(bc, y[:, 0], y[:, -1],
  274. p, bc0)
  275. else:
  276. dbc_dya, dbc_dyb, dbc_dp = bc_jac(y[:, 0], y[:, -1], p)
  277. return construct_global_jac(n, m, k, i_jac, j_jac, h, df_dy,
  278. df_dy_middle, df_dp, df_dp_middle, dbc_dya,
  279. dbc_dyb, dbc_dp)
  280. return col_fun, sys_jac
  281. def solve_newton(n, m, h, col_fun, bc, jac, y, p, B, bvp_tol):
  282. """Solve the nonlinear collocation system by a Newton method.
  283. This is a simple Newton method with a backtracking line search. As
  284. advised in [1]_, an affine-invariant criterion function F = ||J^-1 r||^2
  285. is used, where J is the Jacobian matrix at the current iteration and r is
  286. the vector or collocation residuals (values of the system lhs).
  287. The method alters between full Newton iterations and the fixed-Jacobian
  288. iterations based
  289. There are other tricks proposed in [1]_, but they are not used as they
  290. don't seem to improve anything significantly, and even break the
  291. convergence on some test problems I tried.
  292. All important parameters of the algorithm are defined inside the function.
  293. Parameters
  294. ----------
  295. n : int
  296. Number of equations in the ODE system.
  297. m : int
  298. Number of nodes in the mesh.
  299. h : ndarray, shape (m-1,)
  300. Mesh intervals.
  301. col_fun : callable
  302. Function computing collocation residuals.
  303. bc : callable
  304. Function computing boundary condition residuals.
  305. jac : callable
  306. Function computing the Jacobian of the whole system (including
  307. collocation and boundary condition residuals). It is supposed to
  308. return csc_matrix.
  309. y : ndarray, shape (n, m)
  310. Initial guess for the function values at the mesh nodes.
  311. p : ndarray, shape (k,)
  312. Initial guess for the unknown parameters.
  313. B : ndarray with shape (n, n) or None
  314. Matrix to force the S y(a) = 0 condition for a problems with the
  315. singular term. If None, the singular term is assumed to be absent.
  316. bvp_tol : float
  317. Tolerance to which we want to solve a BVP.
  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. # We also need to control residuals of the boundary conditions. But it
  342. # seems that they become very small eventually as the solver progresses,
  343. # i. e. the tolerance for BC are not very important. We set it 1.5 orders
  344. # lower than the BVP tolerance as well.
  345. tol_bc = 5e-2 * bvp_tol
  346. # Maximum allowed number of Jacobian evaluation and factorization, in
  347. # other words the maximum number of full Newton iterations. A small value
  348. # is recommended in the literature.
  349. max_njev = 4
  350. # Maximum number of iterations, considering that some of them can be
  351. # performed with the fixed Jacobian. In theory such iterations are cheap,
  352. # but it's not that simple in Python.
  353. max_iter = 8
  354. # Minimum relative improvement of the criterion function to accept the
  355. # step (Armijo constant).
  356. sigma = 0.2
  357. # Step size decrease factor for backtracking.
  358. tau = 0.5
  359. # Maximum number of backtracking steps, the minimum step is then
  360. # tau ** n_trial.
  361. n_trial = 4
  362. col_res, y_middle, f, f_middle = col_fun(y, p)
  363. bc_res = bc(y[:, 0], y[:, -1], p)
  364. res = np.hstack((col_res.ravel(order='F'), bc_res))
  365. njev = 0
  366. singular = False
  367. recompute_jac = True
  368. for iteration in range(max_iter):
  369. if recompute_jac:
  370. J = jac(y, p, y_middle, f, f_middle, bc_res)
  371. njev += 1
  372. try:
  373. LU = splu(J)
  374. except RuntimeError:
  375. singular = True
  376. break
  377. step = LU.solve(res)
  378. cost = np.dot(step, step)
  379. y_step = step[:m * n].reshape((n, m), order='F')
  380. p_step = step[m * n:]
  381. alpha = 1
  382. for trial in range(n_trial + 1):
  383. y_new = y - alpha * y_step
  384. if B is not None:
  385. y_new[:, 0] = np.dot(B, y_new[:, 0])
  386. p_new = p - alpha * p_step
  387. col_res, y_middle, f, f_middle = col_fun(y_new, p_new)
  388. bc_res = bc(y_new[:, 0], y_new[:, -1], p_new)
  389. res = np.hstack((col_res.ravel(order='F'), bc_res))
  390. step_new = LU.solve(res)
  391. cost_new = np.dot(step_new, step_new)
  392. if cost_new < (1 - 2 * alpha * sigma) * cost:
  393. break
  394. if trial < n_trial:
  395. alpha *= tau
  396. y = y_new
  397. p = p_new
  398. if njev == max_njev:
  399. break
  400. if (np.all(np.abs(col_res) < tol_r * (1 + np.abs(f_middle))) and
  401. np.all(bc_res < tol_bc)):
  402. break
  403. # If the full step was taken, then we are going to continue with
  404. # the same Jacobian. This is the approach of BVP_SOLVER.
  405. if alpha == 1:
  406. step = step_new
  407. cost = cost_new
  408. recompute_jac = False
  409. else:
  410. recompute_jac = True
  411. return y, p, singular
  412. def print_iteration_header():
  413. print("{:^15}{:^15}{:^15}{:^15}".format(
  414. "Iteration", "Max residual", "Total nodes", "Nodes added"))
  415. def print_iteration_progress(iteration, residual, total_nodes, nodes_added):
  416. print("{:^15}{:^15.2e}{:^15}{:^15}".format(
  417. iteration, residual, total_nodes, nodes_added))
  418. class BVPResult(OptimizeResult):
  419. pass
  420. TERMINATION_MESSAGES = {
  421. 0: "The algorithm converged to the desired accuracy.",
  422. 1: "The maximum number of mesh nodes is exceeded.",
  423. 2: "A singular Jacobian encountered when solving the collocation system."
  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.rollaxis(c, 1)
  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. def solve_bvp(fun, bc, x, y, p=None, S=None, fun_jac=None, bc_jac=None,
  574. tol=1e-3, max_nodes=1000, verbose=0):
  575. """Solve a boundary-value problem for a system of ODEs.
  576. This function numerically solves a first order system of ODEs subject to
  577. two-point boundary conditions::
  578. dy / dx = f(x, y, p) + S * y / (x - a), a <= x <= b
  579. bc(y(a), y(b), p) = 0
  580. Here x is a 1-dimensional independent variable, y(x) is a n-dimensional
  581. vector-valued function and p is a k-dimensional vector of unknown
  582. parameters which is to be found along with y(x). For the problem to be
  583. determined there must be n + k boundary conditions, i.e. bc must be
  584. (n + k)-dimensional function.
  585. The last singular term in the right-hand side of the system is optional.
  586. It is defined by an n-by-n matrix S, such that the solution must satisfy
  587. S y(a) = 0. This condition will be forced during iterations, so it must not
  588. contradict boundary conditions. See [2]_ for the explanation how this term
  589. is handled when solving BVPs numerically.
  590. Problems in a complex domain can be solved as well. In this case y and p
  591. are considered to be complex, and f and bc are assumed to be complex-valued
  592. functions, but x stays real. Note that f and bc must be complex
  593. differentiable (satisfy Cauchy-Riemann equations [4]_), otherwise you
  594. should rewrite your problem for real and imaginary parts separately. To
  595. solve a problem in a complex domain, pass an initial guess for y with a
  596. complex data type (see below).
  597. Parameters
  598. ----------
  599. fun : callable
  600. Right-hand side of the system. The calling signature is ``fun(x, y)``,
  601. or ``fun(x, y, p)`` if parameters are present. All arguments are
  602. ndarray: ``x`` with shape (m,), ``y`` with shape (n, m), meaning that
  603. ``y[:, i]`` corresponds to ``x[i]``, and ``p`` with shape (k,). The
  604. return value must be an array with shape (n, m) and with the same
  605. layout as ``y``.
  606. bc : callable
  607. Function evaluating residuals of the boundary conditions. The calling
  608. signature is ``bc(ya, yb)``, or ``bc(ya, yb, p)`` if parameters are
  609. present. All arguments are ndarray: ``ya`` and ``yb`` with shape (n,),
  610. and ``p`` with shape (k,). The return value must be an array with
  611. shape (n + k,).
  612. x : array_like, shape (m,)
  613. Initial mesh. Must be a strictly increasing sequence of real numbers
  614. with ``x[0]=a`` and ``x[-1]=b``.
  615. y : array_like, shape (n, m)
  616. Initial guess for the function values at the mesh nodes, i-th column
  617. corresponds to ``x[i]``. For problems in a complex domain pass `y`
  618. with a complex data type (even if the initial guess is purely real).
  619. p : array_like with shape (k,) or None, optional
  620. Initial guess for the unknown parameters. If None (default), it is
  621. assumed that the problem doesn't depend on any parameters.
  622. S : array_like with shape (n, n) or None
  623. Matrix defining the singular term. If None (default), the problem is
  624. solved without the singular term.
  625. fun_jac : callable or None, optional
  626. Function computing derivatives of f with respect to y and p. The
  627. calling signature is ``fun_jac(x, y)``, or ``fun_jac(x, y, p)`` if
  628. parameters are present. The return must contain 1 or 2 elements in the
  629. following order:
  630. * df_dy : array_like with shape (n, n, m) where an element
  631. (i, j, q) equals to d f_i(x_q, y_q, p) / d (y_q)_j.
  632. * df_dp : array_like with shape (n, k, m) where an element
  633. (i, j, q) equals to d f_i(x_q, y_q, p) / d p_j.
  634. Here q numbers nodes at which x and y are defined, whereas i and j
  635. number vector components. If the problem is solved without unknown
  636. parameters df_dp should not be returned.
  637. If `fun_jac` is None (default), the derivatives will be estimated
  638. by the forward finite differences.
  639. bc_jac : callable or None, optional
  640. Function computing derivatives of bc with respect to ya, yb and p.
  641. The calling signature is ``bc_jac(ya, yb)``, or ``bc_jac(ya, yb, p)``
  642. if parameters are present. The return must contain 2 or 3 elements in
  643. the following order:
  644. * dbc_dya : array_like with shape (n, n) where an element (i, j)
  645. equals to d bc_i(ya, yb, p) / d ya_j.
  646. * dbc_dyb : array_like with shape (n, n) where an element (i, j)
  647. equals to d bc_i(ya, yb, p) / d yb_j.
  648. * dbc_dp : array_like with shape (n, k) where an element (i, j)
  649. equals to d bc_i(ya, yb, p) / d p_j.
  650. If the problem is solved without unknown parameters dbc_dp should not
  651. be returned.
  652. If `bc_jac` is None (default), the derivatives will be estimated by
  653. the forward finite differences.
  654. tol : float, optional
  655. Desired tolerance of the solution. If we define ``r = y' - f(x, y)``
  656. where y is the found solution, then the solver tries to achieve on each
  657. mesh interval ``norm(r / (1 + abs(f)) < tol``, where ``norm`` is
  658. estimated in a root mean squared sense (using a numerical quadrature
  659. formula). Default is 1e-3.
  660. max_nodes : int, optional
  661. Maximum allowed number of the mesh nodes. If exceeded, the algorithm
  662. terminates. Default is 1000.
  663. verbose : {0, 1, 2}, optional
  664. Level of algorithm's verbosity:
  665. * 0 (default) : work silently.
  666. * 1 : display a termination report.
  667. * 2 : display progress during iterations.
  668. Returns
  669. -------
  670. Bunch object with the following fields defined:
  671. sol : PPoly
  672. Found solution for y as `scipy.interpolate.PPoly` instance, a C1
  673. continuous cubic spline.
  674. p : ndarray or None, shape (k,)
  675. Found parameters. None, if the parameters were not present in the
  676. problem.
  677. x : ndarray, shape (m,)
  678. Nodes of the final mesh.
  679. y : ndarray, shape (n, m)
  680. Solution values at the mesh nodes.
  681. yp : ndarray, shape (n, m)
  682. Solution derivatives at the mesh nodes.
  683. rms_residuals : ndarray, shape (m - 1,)
  684. RMS values of the relative residuals over each mesh interval (see the
  685. description of `tol` parameter).
  686. niter : int
  687. Number of completed iterations.
  688. status : int
  689. Reason for algorithm termination:
  690. * 0: The algorithm converged to the desired accuracy.
  691. * 1: The maximum number of mesh nodes is exceeded.
  692. * 2: A singular Jacobian encountered when solving the collocation
  693. system.
  694. message : string
  695. Verbal description of the termination reason.
  696. success : bool
  697. True if the algorithm converged to the desired accuracy (``status=0``).
  698. Notes
  699. -----
  700. This function implements a 4-th order collocation algorithm with the
  701. control of residuals similar to [1]_. A collocation system is solved
  702. by a damped Newton method with an affine-invariant criterion function as
  703. described in [3]_.
  704. Note that in [1]_ integral residuals are defined without normalization
  705. by interval lengths. So their definition is different by a multiplier of
  706. h**0.5 (h is an interval length) from the definition used here.
  707. .. versionadded:: 0.18.0
  708. References
  709. ----------
  710. .. [1] J. Kierzenka, L. F. Shampine, "A BVP Solver Based on Residual
  711. Control and the Maltab PSE", ACM Trans. Math. Softw., Vol. 27,
  712. Number 3, pp. 299-316, 2001.
  713. .. [2] L.F. Shampine, P. H. Muir and H. Xu, "A User-Friendly Fortran BVP
  714. Solver".
  715. .. [3] U. Ascher, R. Mattheij and R. Russell "Numerical Solution of
  716. Boundary Value Problems for Ordinary Differential Equations".
  717. .. [4] `Cauchy-Riemann equations
  718. <https://en.wikipedia.org/wiki/Cauchy-Riemann_equations>`_ on
  719. Wikipedia.
  720. Examples
  721. --------
  722. In the first example we solve Bratu's problem::
  723. y'' + k * exp(y) = 0
  724. y(0) = y(1) = 0
  725. for k = 1.
  726. We rewrite the equation as a first order system and implement its
  727. right-hand side evaluation::
  728. y1' = y2
  729. y2' = -exp(y1)
  730. >>> def fun(x, y):
  731. ... return np.vstack((y[1], -np.exp(y[0])))
  732. Implement evaluation of the boundary condition residuals:
  733. >>> def bc(ya, yb):
  734. ... return np.array([ya[0], yb[0]])
  735. Define the initial mesh with 5 nodes:
  736. >>> x = np.linspace(0, 1, 5)
  737. This problem is known to have two solutions. To obtain both of them we
  738. use two different initial guesses for y. We denote them by subscripts
  739. a and b.
  740. >>> y_a = np.zeros((2, x.size))
  741. >>> y_b = np.zeros((2, x.size))
  742. >>> y_b[0] = 3
  743. Now we are ready to run the solver.
  744. >>> from scipy.integrate import solve_bvp
  745. >>> res_a = solve_bvp(fun, bc, x, y_a)
  746. >>> res_b = solve_bvp(fun, bc, x, y_b)
  747. Let's plot the two found solutions. We take an advantage of having the
  748. solution in a spline form to produce a smooth plot.
  749. >>> x_plot = np.linspace(0, 1, 100)
  750. >>> y_plot_a = res_a.sol(x_plot)[0]
  751. >>> y_plot_b = res_b.sol(x_plot)[0]
  752. >>> import matplotlib.pyplot as plt
  753. >>> plt.plot(x_plot, y_plot_a, label='y_a')
  754. >>> plt.plot(x_plot, y_plot_b, label='y_b')
  755. >>> plt.legend()
  756. >>> plt.xlabel("x")
  757. >>> plt.ylabel("y")
  758. >>> plt.show()
  759. We see that the two solutions have similar shape, but differ in scale
  760. significantly.
  761. In the second example we solve a simple Sturm-Liouville problem::
  762. y'' + k**2 * y = 0
  763. y(0) = y(1) = 0
  764. It is known that a non-trivial solution y = A * sin(k * x) is possible for
  765. k = pi * n, where n is an integer. To establish the normalization constant
  766. A = 1 we add a boundary condition::
  767. y'(0) = k
  768. Again we rewrite our equation as a first order system and implement its
  769. right-hand side evaluation::
  770. y1' = y2
  771. y2' = -k**2 * y1
  772. >>> def fun(x, y, p):
  773. ... k = p[0]
  774. ... return np.vstack((y[1], -k**2 * y[0]))
  775. Note that parameters p are passed as a vector (with one element in our
  776. case).
  777. Implement the boundary conditions:
  778. >>> def bc(ya, yb, p):
  779. ... k = p[0]
  780. ... return np.array([ya[0], yb[0], ya[1] - k])
  781. Setup the initial mesh and guess for y. We aim to find the solution for
  782. k = 2 * pi, to achieve that we set values of y to approximately follow
  783. sin(2 * pi * x):
  784. >>> x = np.linspace(0, 1, 5)
  785. >>> y = np.zeros((2, x.size))
  786. >>> y[0, 1] = 1
  787. >>> y[0, 3] = -1
  788. Run the solver with 6 as an initial guess for k.
  789. >>> sol = solve_bvp(fun, bc, x, y, p=[6])
  790. We see that the found k is approximately correct:
  791. >>> sol.p[0]
  792. 6.28329460046
  793. And finally plot the solution to see the anticipated sinusoid:
  794. >>> x_plot = np.linspace(0, 1, 100)
  795. >>> y_plot = sol.sol(x_plot)[0]
  796. >>> plt.plot(x_plot, y_plot)
  797. >>> plt.xlabel("x")
  798. >>> plt.ylabel("y")
  799. >>> plt.show()
  800. """
  801. x = np.asarray(x, dtype=float)
  802. if x.ndim != 1:
  803. raise ValueError("`x` must be 1 dimensional.")
  804. h = np.diff(x)
  805. if np.any(h <= 0):
  806. raise ValueError("`x` must be strictly increasing.")
  807. a = x[0]
  808. y = np.asarray(y)
  809. if np.issubdtype(y.dtype, np.complexfloating):
  810. dtype = complex
  811. else:
  812. dtype = float
  813. y = y.astype(dtype, copy=False)
  814. if y.ndim != 2:
  815. raise ValueError("`y` must be 2 dimensional.")
  816. if y.shape[1] != x.shape[0]:
  817. raise ValueError("`y` is expected to have {} columns, but actually "
  818. "has {}.".format(x.shape[0], y.shape[1]))
  819. if p is None:
  820. p = np.array([])
  821. else:
  822. p = np.asarray(p, dtype=dtype)
  823. if p.ndim != 1:
  824. raise ValueError("`p` must be 1 dimensional.")
  825. if tol < 100 * EPS:
  826. warn("`tol` is too low, setting to {:.2e}".format(100 * EPS))
  827. tol = 100 * EPS
  828. if verbose not in [0, 1, 2]:
  829. raise ValueError("`verbose` must be in [0, 1, 2].")
  830. n = y.shape[0]
  831. k = p.shape[0]
  832. if S is not None:
  833. S = np.asarray(S, dtype=dtype)
  834. if S.shape != (n, n):
  835. raise ValueError("`S` is expected to have shape {}, "
  836. "but actually has {}".format((n, n), S.shape))
  837. # Compute I - S^+ S to impose necessary boundary conditions.
  838. B = np.identity(n) - np.dot(pinv(S), S)
  839. y[:, 0] = np.dot(B, y[:, 0])
  840. # Compute (I - S)^+ to correct derivatives at x=a.
  841. D = pinv(np.identity(n) - S)
  842. else:
  843. B = None
  844. D = None
  845. fun_wrapped, bc_wrapped, fun_jac_wrapped, bc_jac_wrapped = wrap_functions(
  846. fun, bc, fun_jac, bc_jac, k, a, S, D, dtype)
  847. f = fun_wrapped(x, y, p)
  848. if f.shape != y.shape:
  849. raise ValueError("`fun` return is expected to have shape {}, "
  850. "but actually has {}.".format(y.shape, f.shape))
  851. bc_res = bc_wrapped(y[:, 0], y[:, -1], p)
  852. if bc_res.shape != (n + k,):
  853. raise ValueError("`bc` return is expected to have shape {}, "
  854. "but actually has {}.".format((n + k,), bc_res.shape))
  855. status = 0
  856. iteration = 0
  857. if verbose == 2:
  858. print_iteration_header()
  859. while True:
  860. m = x.shape[0]
  861. col_fun, jac_sys = prepare_sys(n, m, k, fun_wrapped, bc_wrapped,
  862. fun_jac_wrapped, bc_jac_wrapped, x, h)
  863. y, p, singular = solve_newton(n, m, h, col_fun, bc_wrapped, jac_sys,
  864. y, p, B, tol)
  865. iteration += 1
  866. col_res, y_middle, f, f_middle = collocation_fun(fun_wrapped, y,
  867. p, x, h)
  868. # This relation is not trivial, but can be verified.
  869. r_middle = 1.5 * col_res / h
  870. sol = create_spline(y, f, x, h)
  871. rms_res = estimate_rms_residuals(fun_wrapped, sol, x, h, p,
  872. r_middle, f_middle)
  873. max_rms_res = np.max(rms_res)
  874. if singular:
  875. status = 2
  876. break
  877. insert_1, = np.nonzero((rms_res > tol) & (rms_res < 100 * tol))
  878. insert_2, = np.nonzero(rms_res >= 100 * tol)
  879. nodes_added = insert_1.shape[0] + 2 * insert_2.shape[0]
  880. if m + nodes_added > max_nodes:
  881. status = 1
  882. if verbose == 2:
  883. nodes_added = "({})".format(nodes_added)
  884. print_iteration_progress(iteration, max_rms_res, m,
  885. nodes_added)
  886. break
  887. if verbose == 2:
  888. print_iteration_progress(iteration, max_rms_res, m, nodes_added)
  889. if nodes_added > 0:
  890. x = modify_mesh(x, insert_1, insert_2)
  891. h = np.diff(x)
  892. y = sol(x)
  893. else:
  894. status = 0
  895. break
  896. if verbose > 0:
  897. if status == 0:
  898. print("Solved in {} iterations, number of nodes {}, "
  899. "maximum relative residual {:.2e}."
  900. .format(iteration, x.shape[0], max_rms_res))
  901. elif status == 1:
  902. print("Number of nodes is exceeded after iteration {}, "
  903. "maximum relative residual {:.2e}."
  904. .format(iteration, max_rms_res))
  905. elif status == 2:
  906. print("Singular Jacobian encountered when solving the collocation "
  907. "system on iteration {}, maximum relative residual {:.2e}."
  908. .format(iteration, max_rms_res))
  909. if p.size == 0:
  910. p = None
  911. return BVPResult(sol=sol, p=p, x=x, y=y, yp=f, rms_residuals=rms_res,
  912. niter=iteration, status=status,
  913. message=TERMINATION_MESSAGES[status], success=status == 0)