gsl_ode-initval__gear1.c 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273
  1. /* ode-initval/gear1.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. /* Gear 1. This is the implicit Euler a.k.a backward Euler method. */
  20. /* Author: G. Jungman
  21. */
  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 *k;
  37. double *y0;
  38. double *y0_orig;
  39. double *y_onestep;
  40. }
  41. gear1_state_t;
  42. static void *
  43. gear1_alloc (size_t dim)
  44. {
  45. gear1_state_t *state = (gear1_state_t *) malloc (sizeof (gear1_state_t));
  46. if (state == 0)
  47. {
  48. GSL_ERROR_NULL ("failed to allocate space for gear1_state", GSL_ENOMEM);
  49. }
  50. state->k = (double *) malloc (dim * sizeof (double));
  51. if (state->k == 0)
  52. {
  53. free (state);
  54. GSL_ERROR_NULL ("failed to allocate space for k", GSL_ENOMEM);
  55. }
  56. state->y0 = (double *) malloc (dim * sizeof (double));
  57. if (state->y0 == 0)
  58. {
  59. free (state->k);
  60. free (state);
  61. GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
  62. }
  63. state->y0_orig = (double *) malloc (dim * sizeof (double));
  64. if (state->y0_orig == 0)
  65. {
  66. free (state->y0);
  67. free (state->k);
  68. free (state);
  69. GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
  70. }
  71. state->y_onestep = (double *) malloc (dim * sizeof (double));
  72. if (state->y_onestep == 0)
  73. {
  74. free (state->y0_orig);
  75. free (state->y0);
  76. free (state->k);
  77. free (state);
  78. GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
  79. }
  80. return state;
  81. }
  82. static int
  83. gear1_step (double *y, gear1_state_t *state,
  84. const double h, const double t,
  85. const size_t dim, const gsl_odeiv_system *sys)
  86. {
  87. /* Makes an implicit Euler advance with step size h.
  88. y0 is the initial values of variables y.
  89. The implicit matrix equations to solve are:
  90. k = y0 + h * f(t + h, k)
  91. y = y0 + h * f(t + h, k)
  92. */
  93. const int iter_steps = 3;
  94. int nu;
  95. size_t i;
  96. double *y0 = state->y0;
  97. double *k = state->k;
  98. /* Iterative solution of k = y0 + h * f(t + h, k)
  99. Note: This method does not check for convergence of the
  100. iterative solution!
  101. */
  102. for (nu = 0; nu < iter_steps; nu++)
  103. {
  104. int s = GSL_ODEIV_FN_EVAL(sys, t + h, y, k);
  105. if (s != GSL_SUCCESS)
  106. {
  107. return s;
  108. }
  109. for (i=0; i<dim; i++)
  110. {
  111. y[i] = y0[i] + h * k[i];
  112. }
  113. }
  114. return GSL_SUCCESS;
  115. }
  116. static int
  117. gear1_apply(void * vstate,
  118. size_t dim,
  119. double t,
  120. double h,
  121. double y[],
  122. double yerr[],
  123. const double dydt_in[],
  124. double dydt_out[],
  125. const gsl_odeiv_system * sys)
  126. {
  127. gear1_state_t *state = (gear1_state_t *) vstate;
  128. size_t i;
  129. double *y0 = state->y0;
  130. double *y0_orig = state->y0_orig;
  131. double *y_onestep = state->y_onestep;
  132. /* initialization */
  133. DBL_MEMCPY(y0, y, dim);
  134. /* Save initial values for possible failures */
  135. DBL_MEMCPY (y0_orig, y, dim);
  136. /* First traverse h with one step (save to y_onestep) */
  137. DBL_MEMCPY (y_onestep, y, dim);
  138. {
  139. int s = gear1_step (y_onestep, state, h, t, dim, sys);
  140. if (s != GSL_SUCCESS)
  141. {
  142. return s;
  143. }
  144. }
  145. /* Then with two steps with half step length (save to y) */
  146. {
  147. int s = gear1_step (y, state, h / 2.0, t, dim, sys);
  148. if (s != GSL_SUCCESS)
  149. {
  150. /* Restore original y vector */
  151. DBL_MEMCPY (y, y0_orig, dim);
  152. return s;
  153. }
  154. }
  155. DBL_MEMCPY (y0, y, dim);
  156. {
  157. int s = gear1_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
  158. if (s != GSL_SUCCESS)
  159. {
  160. /* Restore original y vector */
  161. DBL_MEMCPY (y, y0_orig, dim);
  162. return s;
  163. }
  164. }
  165. /* Cleanup update */
  166. if (dydt_out != NULL)
  167. {
  168. int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
  169. if (s != GSL_SUCCESS)
  170. {
  171. /* Restore original y vector */
  172. DBL_MEMCPY (y, y0_orig, dim);
  173. return s;
  174. }
  175. }
  176. /* Error estimation */
  177. for (i = 0; i < dim; i++)
  178. {
  179. yerr[i] = 4.0 * (y[i] - y_onestep[i]);
  180. }
  181. return GSL_SUCCESS;
  182. }
  183. static int
  184. gear1_reset (void *vstate, size_t dim)
  185. {
  186. gear1_state_t *state = (gear1_state_t *) vstate;
  187. DBL_ZERO_MEMSET (state->y_onestep, dim);
  188. DBL_ZERO_MEMSET (state->y0_orig, dim);
  189. DBL_ZERO_MEMSET (state->y0, dim);
  190. DBL_ZERO_MEMSET (state->k, dim);
  191. return GSL_SUCCESS;
  192. }
  193. static unsigned int
  194. gear1_order (void *vstate)
  195. {
  196. gear1_state_t *state = (gear1_state_t *) vstate;
  197. state = 0; /* prevent warnings about unused parameters */
  198. return 1;
  199. }
  200. static void
  201. gear1_free (void *vstate)
  202. {
  203. gear1_state_t *state = (gear1_state_t *) vstate;
  204. free (state->y_onestep);
  205. free (state->y0_orig);
  206. free (state->y0);
  207. free (state->k);
  208. free (state);
  209. }
  210. static const gsl_odeiv_step_type gear1_type = { "gear1", /* name */
  211. 0, /* can use dydt_in? */
  212. 1, /* gives exact dydt_out? */
  213. &gear1_alloc,
  214. &gear1_apply,
  215. &gear1_reset,
  216. &gear1_order,
  217. &gear1_free
  218. };
  219. const gsl_odeiv_step_type *gsl_odeiv_step_gear1 = &gear1_type;