gsl_ode-initval__rk2imp.c 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337
  1. /* ode-initval/rk2imp.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
  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 the implicit
  20. midpoint rule. */
  21. /* Author: G. Jungman */
  22. /* Error estimation by step doubling, see eg. Ascher, U.M., Petzold,
  23. L.R., Computer methods for ordinary differential and
  24. differential-algebraic equations, SIAM, Philadelphia, 1998.
  25. The method is also described in eg. this reference.
  26. */
  27. #include "gsl__config.h"
  28. #include <stdlib.h>
  29. #include <string.h>
  30. #include "gsl_math.h"
  31. #include "gsl_errno.h"
  32. #include "gsl_odeiv.h"
  33. #include "gsl_ode-initval__odeiv_util.h"
  34. typedef struct
  35. {
  36. double *Y1;
  37. double *y0;
  38. double *ytmp;
  39. double *y_onestep;
  40. double *y0_orig;
  41. }
  42. rk2imp_state_t;
  43. static void *
  44. rk2imp_alloc (size_t dim)
  45. {
  46. rk2imp_state_t *state = (rk2imp_state_t *) malloc (sizeof (rk2imp_state_t));
  47. if (state == 0)
  48. {
  49. GSL_ERROR_NULL ("failed to allocate space for rk2imp_state",
  50. GSL_ENOMEM);
  51. }
  52. state->Y1 = (double *) malloc (dim * sizeof (double));
  53. if (state->Y1 == 0)
  54. {
  55. free (state);
  56. GSL_ERROR_NULL ("failed to allocate space for Y1", GSL_ENOMEM);
  57. }
  58. state->ytmp = (double *) malloc (dim * sizeof (double));
  59. if (state->ytmp == 0)
  60. {
  61. free (state->Y1);
  62. free (state);
  63. GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
  64. }
  65. state->y0 = (double *) malloc (dim * sizeof (double));
  66. if (state->y0 == 0)
  67. {
  68. free (state->Y1);
  69. free (state->ytmp);
  70. free (state);
  71. GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
  72. }
  73. state->y_onestep = (double *) malloc (dim * sizeof (double));
  74. if (state->y_onestep == 0)
  75. {
  76. free (state->Y1);
  77. free (state->ytmp);
  78. free (state->y0);
  79. free (state);
  80. GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
  81. }
  82. state->y0_orig = (double *) malloc (dim * sizeof (double));
  83. if (state->y0_orig == 0)
  84. {
  85. free (state->y_onestep);
  86. free (state->Y1);
  87. free (state->ytmp);
  88. free (state->y0);
  89. free (state);
  90. GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
  91. }
  92. return state;
  93. }
  94. static int
  95. rk2imp_step (double *y, rk2imp_state_t *state,
  96. const double h, const double t,
  97. const size_t dim, const gsl_odeiv_system *sys)
  98. {
  99. /* Makes a Runge-Kutta 2nd order implicit advance with step size h.
  100. y0 is initial values of variables y.
  101. The implicit matrix equations to solve are:
  102. Y1 = y0 + h/2 * f(t + h/2, Y1)
  103. y = y0 + h * f(t + h/2, Y1)
  104. */
  105. const double *y0 = state->y0;
  106. double *Y1 = state->Y1;
  107. double *ytmp = state->ytmp;
  108. int max_iter=3;
  109. int nu;
  110. size_t i;
  111. /* iterative solution of Y1 = y0 + h/2 * f(t + h/2, Y1)
  112. Y1 should include initial values at call.
  113. Note: This method does not check for convergence of the
  114. iterative solution!
  115. */
  116. for (nu = 0; nu < max_iter; nu++)
  117. {
  118. for (i = 0; i < dim; i++)
  119. {
  120. ytmp[i] = y0[i] + 0.5 * h * Y1[i];
  121. }
  122. {
  123. int s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, ytmp, Y1);
  124. if (s != GSL_SUCCESS)
  125. {
  126. return s;
  127. }
  128. }
  129. }
  130. /* assignment */
  131. for (i = 0; i < dim; i++)
  132. {
  133. y[i] = y0[i] + h * Y1[i];
  134. }
  135. return GSL_SUCCESS;
  136. }
  137. static int
  138. rk2imp_apply (void *vstate,
  139. size_t dim,
  140. double t,
  141. double h,
  142. double y[],
  143. double yerr[],
  144. const double dydt_in[],
  145. double dydt_out[], const gsl_odeiv_system * sys)
  146. {
  147. rk2imp_state_t *state = (rk2imp_state_t *) vstate;
  148. size_t i;
  149. double *Y1 = state->Y1;
  150. double *y0 = state->y0;
  151. double *y_onestep = state->y_onestep;
  152. double *y0_orig = state->y0_orig;
  153. /* Error estimation is done by step doubling procedure */
  154. /* initialization step */
  155. DBL_MEMCPY (y0, y, dim);
  156. /* Save initial values for possible failures */
  157. DBL_MEMCPY (y0_orig, y, dim);
  158. if (dydt_in != NULL)
  159. {
  160. DBL_MEMCPY (Y1, dydt_in, dim);
  161. }
  162. else
  163. {
  164. int s = GSL_ODEIV_FN_EVAL (sys, t, y, Y1);
  165. if (s != GSL_SUCCESS)
  166. {
  167. return s;
  168. }
  169. }
  170. /* First traverse h with one step (save to y_onestep) */
  171. DBL_MEMCPY (y_onestep, y, dim);
  172. {
  173. int s = rk2imp_step (y_onestep, state, h, t, dim, sys);
  174. if (s != GSL_SUCCESS)
  175. {
  176. return s;
  177. }
  178. }
  179. /* Then with two steps with half step length (save to y) */
  180. {
  181. int s = rk2imp_step (y, state, h / 2.0, t, dim, sys);
  182. if (s != GSL_SUCCESS)
  183. {
  184. /* Restore original y vector */
  185. DBL_MEMCPY (y, y0_orig, dim);
  186. return s;
  187. }
  188. }
  189. DBL_MEMCPY (y0, y, dim);
  190. {
  191. int s = GSL_ODEIV_FN_EVAL (sys, t + h / 2.0, y, Y1);
  192. if (s != GSL_SUCCESS)
  193. {
  194. /* Restore original y vector */
  195. DBL_MEMCPY (y, y0_orig, dim);
  196. return s;
  197. }
  198. }
  199. {
  200. int s = rk2imp_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
  201. if (s != GSL_SUCCESS)
  202. {
  203. /* Restore original y vector */
  204. DBL_MEMCPY (y, y0_orig, dim);
  205. return s;
  206. }
  207. }
  208. /* Derivatives at output */
  209. if (dydt_out != NULL)
  210. {
  211. int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
  212. if (s != GSL_SUCCESS)
  213. {
  214. /* Restore original y vector */
  215. DBL_MEMCPY (y, y0_orig, dim);
  216. return s;
  217. }
  218. }
  219. /* Error estimation */
  220. for (i = 0; i < dim; i++)
  221. {
  222. yerr[i] = 4.0 * (y[i] - y_onestep[i]) / 3.0;
  223. }
  224. return GSL_SUCCESS;
  225. }
  226. static int
  227. rk2imp_reset (void *vstate, size_t dim)
  228. {
  229. rk2imp_state_t *state = (rk2imp_state_t *) vstate;
  230. DBL_ZERO_MEMSET (state->Y1, dim);
  231. DBL_ZERO_MEMSET (state->ytmp, dim);
  232. DBL_ZERO_MEMSET (state->y0, dim);
  233. DBL_ZERO_MEMSET (state->y_onestep, dim);
  234. DBL_ZERO_MEMSET (state->y0_orig, dim);
  235. return GSL_SUCCESS;
  236. }
  237. static unsigned int
  238. rk2imp_order (void *vstate)
  239. {
  240. rk2imp_state_t *state = (rk2imp_state_t *) vstate;
  241. state = 0; /* prevent warnings about unused parameters */
  242. return 2;
  243. }
  244. static void
  245. rk2imp_free (void *vstate)
  246. {
  247. rk2imp_state_t *state = (rk2imp_state_t *) vstate;
  248. free (state->Y1);
  249. free (state->ytmp);
  250. free (state->y0);
  251. free (state->y_onestep);
  252. free (state->y0_orig);
  253. free (state);
  254. }
  255. static const gsl_odeiv_step_type rk2imp_type = { "rk2imp", /* name */
  256. 1, /* can use dydt_in */
  257. 1, /* gives exact dydt_out */
  258. &rk2imp_alloc,
  259. &rk2imp_apply,
  260. &rk2imp_reset,
  261. &rk2imp_order,
  262. &rk2imp_free
  263. };
  264. const gsl_odeiv_step_type *gsl_odeiv_step_rk2imp = &rk2imp_type;