gsl_ode-initval__rk2simp.c 8.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383
  1. /* ode-initval/rk2simp.c
  2. *
  3. * Copyright (C) 2004 Tuomo Keskitalo
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. /* Runge-Kutta 2, Gaussian implicit. Also known as implicit midpoint rule.
  20. Non-linear equations solved by linearization, LU-decomposition
  21. and matrix inversion. For reference, see eg.
  22. Ascher, U.M., Petzold, L.R., Computer methods for ordinary
  23. differential and differential-algebraic equations, SIAM,
  24. Philadelphia, 1998.
  25. */
  26. #include "gsl__config.h"
  27. #include <stdlib.h>
  28. #include <string.h>
  29. #include "gsl_math.h"
  30. #include "gsl_errno.h"
  31. #include "gsl_odeiv.h"
  32. #include "gsl_linalg.h"
  33. #include "gsl_ode-initval__odeiv_util.h"
  34. typedef struct
  35. {
  36. double *Y1;
  37. double *y0;
  38. double *y0_orig;
  39. double *ytmp;
  40. double *dfdy; /* Jacobian */
  41. double *dfdt; /* time derivatives, not used */
  42. double *y_onestep;
  43. gsl_permutation *p;
  44. }
  45. rk2simp_state_t;
  46. static void *
  47. rk2simp_alloc (size_t dim)
  48. {
  49. rk2simp_state_t *state =
  50. (rk2simp_state_t *) malloc (sizeof (rk2simp_state_t));
  51. if (state == 0)
  52. {
  53. GSL_ERROR_NULL ("failed to allocate space for rk2simp_state",
  54. GSL_ENOMEM);
  55. }
  56. state->Y1 = (double *) malloc (dim * sizeof (double));
  57. if (state->Y1 == 0)
  58. {
  59. free (state);
  60. GSL_ERROR_NULL ("failed to allocate space for Y1", GSL_ENOMEM);
  61. }
  62. state->y0 = (double *) malloc (dim * sizeof (double));
  63. if (state->y0 == 0)
  64. {
  65. free (state->Y1);
  66. free (state);
  67. GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
  68. }
  69. state->y0_orig = (double *) malloc (dim * sizeof (double));
  70. if (state->y0_orig == 0)
  71. {
  72. free (state->Y1);
  73. free (state->y0);
  74. free (state);
  75. GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
  76. }
  77. state->ytmp = (double *) malloc (dim * sizeof (double));
  78. if (state->ytmp == 0)
  79. {
  80. free (state->Y1);
  81. free (state->y0);
  82. free (state->y0_orig);
  83. free (state);
  84. GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
  85. }
  86. state->dfdy = (double *) malloc (dim * dim * sizeof (double));
  87. if (state->dfdy == 0)
  88. {
  89. free (state->Y1);
  90. free (state->y0);
  91. free (state->y0_orig);
  92. free (state->ytmp);
  93. free (state);
  94. GSL_ERROR_NULL ("failed to allocate space for dfdy", GSL_ENOMEM);
  95. }
  96. state->dfdt = (double *) malloc (dim * sizeof (double));
  97. if (state->dfdt == 0)
  98. {
  99. free (state->Y1);
  100. free (state->y0);
  101. free (state->y0_orig);
  102. free (state->ytmp);
  103. free (state->dfdy);
  104. free (state);
  105. GSL_ERROR_NULL ("failed to allocate space for dfdt", GSL_ENOMEM);
  106. }
  107. state->y_onestep = (double *) malloc (dim * sizeof (double));
  108. if (state->y_onestep == 0)
  109. {
  110. free (state->Y1);
  111. free (state->y0);
  112. free (state->y0_orig);
  113. free (state->ytmp);
  114. free (state->dfdy);
  115. free (state->dfdt);
  116. free (state);
  117. GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
  118. }
  119. state->p = gsl_permutation_alloc (dim);
  120. if (state->p == 0)
  121. {
  122. free (state->Y1);
  123. free (state->y0);
  124. free (state->y0_orig);
  125. free (state->ytmp);
  126. free (state->dfdy);
  127. free (state->dfdt);
  128. free (state);
  129. GSL_ERROR_NULL ("failed to allocate space for p", GSL_ENOMEM);
  130. }
  131. return state;
  132. }
  133. static int
  134. rk2simp_step (double *y, rk2simp_state_t * state,
  135. const double h, const double t,
  136. const size_t dim, const gsl_odeiv_system * sys)
  137. {
  138. /* Makes a Runge-Kutta 2nd order semi-implicit advance with step size h.
  139. y0 is initial values of variables y.
  140. The linearized semi-implicit equations to calculate are:
  141. Y1 = y0 + h/2 * (1 - h/2 * df/dy)^(-1) * f(t + h/2, y0)
  142. y = y0 + h * f(t + h/2, Y1)
  143. */
  144. const double *y0 = state->y0;
  145. double *Y1 = state->Y1;
  146. double *ytmp = state->ytmp;
  147. size_t i;
  148. int s, ps;
  149. gsl_matrix_view J = gsl_matrix_view_array (state->dfdy, dim, dim);
  150. /* First solve Y1.
  151. Calculate the inverse matrix (1 - h/2 * df/dy)^-1
  152. */
  153. /* Create matrix to J */
  154. s = GSL_ODEIV_JA_EVAL (sys, t, y0, state->dfdy, state->dfdt);
  155. if (s != GSL_SUCCESS)
  156. {
  157. return s;
  158. }
  159. gsl_matrix_scale (&J.matrix, -h / 2.0);
  160. gsl_matrix_add_diagonal(&J.matrix, 1.0);
  161. /* Invert it by LU-decomposition to invmat */
  162. s += gsl_linalg_LU_decomp (&J.matrix, state->p, &ps);
  163. if (s != GSL_SUCCESS)
  164. {
  165. return GSL_EFAILED;
  166. }
  167. /* Evaluate f(t + h/2, y0) */
  168. s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, y0, ytmp);
  169. if (s != GSL_SUCCESS)
  170. {
  171. return s;
  172. }
  173. /* Calculate Y1 = y0 + h/2 * ((1-h/2 * df/dy)^-1) ytmp */
  174. {
  175. gsl_vector_const_view y0_view = gsl_vector_const_view_array(y0, dim);
  176. gsl_vector_view ytmp_view = gsl_vector_view_array(ytmp, dim);
  177. gsl_vector_view Y1_view = gsl_vector_view_array(Y1, dim);
  178. s = gsl_linalg_LU_solve (&J.matrix, state->p,
  179. &ytmp_view.vector, &Y1_view.vector);
  180. gsl_vector_scale (&Y1_view.vector, 0.5 * h);
  181. gsl_vector_add (&Y1_view.vector, &y0_view.vector);
  182. }
  183. /* And finally evaluation of f(t + h/2, Y1) and calculation of y */
  184. s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, Y1, ytmp);
  185. if (s != GSL_SUCCESS)
  186. {
  187. return s;
  188. }
  189. for (i = 0; i < dim; i++)
  190. {
  191. y[i] = y0[i] + h * ytmp[i];
  192. }
  193. return s;
  194. }
  195. static int
  196. rk2simp_apply (void *vstate, size_t dim, double t, double h,
  197. double y[], double yerr[], const double dydt_in[],
  198. double dydt_out[], const gsl_odeiv_system * sys)
  199. {
  200. rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  201. size_t i;
  202. double *y0 = state->y0;
  203. double *y0_orig = state->y0_orig;
  204. double *y_onestep = state->y_onestep;
  205. /* Error estimation is done by step doubling procedure */
  206. DBL_MEMCPY (y0, y, dim);
  207. /* Save initial values in case of failure */
  208. DBL_MEMCPY (y0_orig, y, dim);
  209. /* First traverse h with one step (save to y_onestep) */
  210. DBL_MEMCPY (y_onestep, y, dim);
  211. {
  212. int s = rk2simp_step (y_onestep, state, h, t, dim, sys);
  213. if (s != GSL_SUCCESS)
  214. {
  215. return s;
  216. }
  217. }
  218. /* Then with two steps with half step length (save to y) */
  219. {
  220. int s = rk2simp_step (y, state, h / 2.0, t, dim, sys);
  221. if (s != GSL_SUCCESS)
  222. {
  223. /* Restore original y vector */
  224. DBL_MEMCPY (y, y0_orig, dim);
  225. return s;
  226. }
  227. }
  228. DBL_MEMCPY (y0, y, dim);
  229. {
  230. int s = rk2simp_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
  231. if (s != GSL_SUCCESS)
  232. {
  233. /* Restore original y vector */
  234. DBL_MEMCPY (y, y0_orig, dim);
  235. return s;
  236. }
  237. }
  238. /* Derivatives at output */
  239. if (dydt_out != NULL)
  240. {
  241. int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
  242. if (s != GSL_SUCCESS)
  243. {
  244. /* Restore original y vector */
  245. DBL_MEMCPY (y, y0_orig, dim);
  246. return s;
  247. }
  248. }
  249. /* Error estimation */
  250. for (i = 0; i < dim; i++)
  251. {
  252. yerr[i] = 4.0 * (y[i] - y_onestep[i]) / 3.0;
  253. }
  254. return GSL_SUCCESS;
  255. }
  256. static int
  257. rk2simp_reset (void *vstate, size_t dim)
  258. {
  259. rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  260. DBL_ZERO_MEMSET (state->Y1, dim);
  261. DBL_ZERO_MEMSET (state->y0, dim);
  262. DBL_ZERO_MEMSET (state->y0_orig, dim);
  263. DBL_ZERO_MEMSET (state->ytmp, dim);
  264. DBL_ZERO_MEMSET (state->dfdt, dim * dim);
  265. DBL_ZERO_MEMSET (state->dfdt, dim);
  266. DBL_ZERO_MEMSET (state->y_onestep, dim);
  267. return GSL_SUCCESS;
  268. }
  269. static unsigned int
  270. rk2simp_order (void *vstate)
  271. {
  272. rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  273. state = 0; /* prevent warnings about unused parameters */
  274. return 2;
  275. }
  276. static void
  277. rk2simp_free (void *vstate)
  278. {
  279. rk2simp_state_t *state = (rk2simp_state_t *) vstate;
  280. free (state->Y1);
  281. free (state->y0);
  282. free (state->y0_orig);
  283. free (state->ytmp);
  284. free (state->dfdy);
  285. free (state->dfdt);
  286. free (state->y_onestep);
  287. gsl_permutation_free (state->p);
  288. free (state);
  289. }
  290. static const gsl_odeiv_step_type rk2simp_type = {
  291. "rk2simp", /* name */
  292. 0, /* can use dydt_in? */
  293. 1, /* gives exact dydt_out? */
  294. &rk2simp_alloc,
  295. &rk2simp_apply,
  296. &rk2simp_reset,
  297. &rk2simp_order,
  298. &rk2simp_free
  299. };
  300. const gsl_odeiv_step_type *gsl_odeiv_step_rk2simp = &rk2simp_type;