gsl_ode-initval__rkf45.c 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373
  1. /* ode-initval/rkf45.c
  2. *
  3. * Copyright (C) 2001, 2004, 2007 Brian Gough
  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-Fehlberg 4(5)*/
  20. /* Reference eg. Hairer, E., Norsett S.P., Wanner, G. Solving ordinary
  21. differential equations I, Nonstiff Problems, 2nd revised edition,
  22. Springer, 2000.
  23. */
  24. #include "gsl__config.h"
  25. #include <stdlib.h>
  26. #include <string.h>
  27. #include "gsl_errno.h"
  28. #include "gsl_odeiv.h"
  29. #include "gsl_ode-initval__odeiv_util.h"
  30. /* Runge-Kutta-Fehlberg coefficients. Zero elements left out */
  31. static const double ah[] = { 1.0/4.0, 3.0/8.0, 12.0/13.0, 1.0, 1.0/2.0 };
  32. static const double b3[] = { 3.0/32.0, 9.0/32.0 };
  33. static const double b4[] = { 1932.0/2197.0, -7200.0/2197.0, 7296.0/2197.0};
  34. static const double b5[] = { 8341.0/4104.0, -32832.0/4104.0, 29440.0/4104.0, -845.0/4104.0};
  35. static const double b6[] = { -6080.0/20520.0, 41040.0/20520.0, -28352.0/20520.0, 9295.0/20520.0, -5643.0/20520.0};
  36. static const double c1 = 902880.0/7618050.0;
  37. static const double c3 = 3953664.0/7618050.0;
  38. static const double c4 = 3855735.0/7618050.0;
  39. static const double c5 = -1371249.0/7618050.0;
  40. static const double c6 = 277020.0/7618050.0;
  41. /* These are the differences of fifth and fourth order coefficients
  42. for error estimation */
  43. static const double ec[] = { 0.0,
  44. 1.0 / 360.0,
  45. 0.0,
  46. -128.0 / 4275.0,
  47. -2197.0 / 75240.0,
  48. 1.0 / 50.0,
  49. 2.0 / 55.0
  50. };
  51. typedef struct
  52. {
  53. double *k1;
  54. double *k2;
  55. double *k3;
  56. double *k4;
  57. double *k5;
  58. double *k6;
  59. double *y0;
  60. double *ytmp;
  61. }
  62. rkf45_state_t;
  63. static void *
  64. rkf45_alloc (size_t dim)
  65. {
  66. rkf45_state_t *state = (rkf45_state_t *) malloc (sizeof (rkf45_state_t));
  67. if (state == 0)
  68. {
  69. GSL_ERROR_NULL ("failed to allocate space for rkf45_state", GSL_ENOMEM);
  70. }
  71. state->k1 = (double *) malloc (dim * sizeof (double));
  72. if (state->k1 == 0)
  73. {
  74. free (state);
  75. GSL_ERROR_NULL ("failed to allocate space for k1", GSL_ENOMEM);
  76. }
  77. state->k2 = (double *) malloc (dim * sizeof (double));
  78. if (state->k2 == 0)
  79. {
  80. free (state->k1);
  81. free (state);
  82. GSL_ERROR_NULL ("failed to allocate space for k2", GSL_ENOMEM);
  83. }
  84. state->k3 = (double *) malloc (dim * sizeof (double));
  85. if (state->k3 == 0)
  86. {
  87. free (state->k2);
  88. free (state->k1);
  89. free (state);
  90. GSL_ERROR_NULL ("failed to allocate space for k3", GSL_ENOMEM);
  91. }
  92. state->k4 = (double *) malloc (dim * sizeof (double));
  93. if (state->k4 == 0)
  94. {
  95. free (state->k3);
  96. free (state->k2);
  97. free (state->k1);
  98. free (state);
  99. GSL_ERROR_NULL ("failed to allocate space for k4", GSL_ENOMEM);
  100. }
  101. state->k5 = (double *) malloc (dim * sizeof (double));
  102. if (state->k5 == 0)
  103. {
  104. free (state->k4);
  105. free (state->k3);
  106. free (state->k2);
  107. free (state->k1);
  108. free (state);
  109. GSL_ERROR_NULL ("failed to allocate space for k5", GSL_ENOMEM);
  110. }
  111. state->k6 = (double *) malloc (dim * sizeof (double));
  112. if (state->k6 == 0)
  113. {
  114. free (state->k5);
  115. free (state->k4);
  116. free (state->k3);
  117. free (state->k2);
  118. free (state->k1);
  119. free (state);
  120. GSL_ERROR_NULL ("failed to allocate space for k6", GSL_ENOMEM);
  121. }
  122. state->y0 = (double *) malloc (dim * sizeof (double));
  123. if (state->y0 == 0)
  124. {
  125. free (state->k6);
  126. free (state->k5);
  127. free (state->k4);
  128. free (state->k3);
  129. free (state->k2);
  130. free (state->k1);
  131. free (state);
  132. GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
  133. }
  134. state->ytmp = (double *) malloc (dim * sizeof (double));
  135. if (state->ytmp == 0)
  136. {
  137. free (state->y0);
  138. free (state->k6);
  139. free (state->k5);
  140. free (state->k4);
  141. free (state->k3);
  142. free (state->k2);
  143. free (state->k1);
  144. free (state);
  145. GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
  146. }
  147. return state;
  148. }
  149. static int
  150. rkf45_apply (void *vstate,
  151. size_t dim,
  152. double t,
  153. double h,
  154. double y[],
  155. double yerr[],
  156. const double dydt_in[],
  157. double dydt_out[], const gsl_odeiv_system * sys)
  158. {
  159. rkf45_state_t *state = (rkf45_state_t *) vstate;
  160. size_t i;
  161. double *const k1 = state->k1;
  162. double *const k2 = state->k2;
  163. double *const k3 = state->k3;
  164. double *const k4 = state->k4;
  165. double *const k5 = state->k5;
  166. double *const k6 = state->k6;
  167. double *const ytmp = state->ytmp;
  168. double *const y0 = state->y0;
  169. DBL_MEMCPY (y0, y, dim);
  170. /* k1 step */
  171. if (dydt_in != NULL)
  172. {
  173. DBL_MEMCPY (k1, dydt_in, dim);
  174. }
  175. else
  176. {
  177. int s = GSL_ODEIV_FN_EVAL (sys, t, y, k1);
  178. if (s != GSL_SUCCESS)
  179. {
  180. return s;
  181. }
  182. }
  183. for (i = 0; i < dim; i++)
  184. ytmp[i] = y[i] + ah[0] * h * k1[i];
  185. /* k2 step */
  186. {
  187. int s = GSL_ODEIV_FN_EVAL (sys, t + ah[0] * h, ytmp, k2);
  188. if (s != GSL_SUCCESS)
  189. {
  190. return s;
  191. }
  192. }
  193. for (i = 0; i < dim; i++)
  194. ytmp[i] = y[i] + h * (b3[0] * k1[i] + b3[1] * k2[i]);
  195. /* k3 step */
  196. {
  197. int s = GSL_ODEIV_FN_EVAL (sys, t + ah[1] * h, ytmp, k3);
  198. if (s != GSL_SUCCESS)
  199. {
  200. return s;
  201. }
  202. }
  203. for (i = 0; i < dim; i++)
  204. ytmp[i] = y[i] + h * (b4[0] * k1[i] + b4[1] * k2[i] + b4[2] * k3[i]);
  205. /* k4 step */
  206. {
  207. int s = GSL_ODEIV_FN_EVAL (sys, t + ah[2] * h, ytmp, k4);
  208. if (s != GSL_SUCCESS)
  209. {
  210. return s;
  211. }
  212. }
  213. for (i = 0; i < dim; i++)
  214. ytmp[i] =
  215. y[i] + h * (b5[0] * k1[i] + b5[1] * k2[i] + b5[2] * k3[i] +
  216. b5[3] * k4[i]);
  217. /* k5 step */
  218. {
  219. int s = GSL_ODEIV_FN_EVAL (sys, t + ah[3] * h, ytmp, k5);
  220. if (s != GSL_SUCCESS)
  221. {
  222. return s;
  223. }
  224. }
  225. for (i = 0; i < dim; i++)
  226. ytmp[i] =
  227. y[i] + h * (b6[0] * k1[i] + b6[1] * k2[i] + b6[2] * k3[i] +
  228. b6[3] * k4[i] + b6[4] * k5[i]);
  229. /* k6 step and final sum */
  230. {
  231. int s = GSL_ODEIV_FN_EVAL (sys, t + ah[4] * h, ytmp, k6);
  232. if (s != GSL_SUCCESS)
  233. {
  234. return s;
  235. }
  236. }
  237. for (i = 0; i < dim; i++)
  238. {
  239. const double d_i = c1 * k1[i] + c3 * k3[i] + c4 * k4[i] + c5 * k5[i] + c6 * k6[i];
  240. y[i] += h * d_i;
  241. }
  242. /* Derivatives at output */
  243. if (dydt_out != NULL)
  244. {
  245. int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
  246. if (s != GSL_SUCCESS)
  247. {
  248. /* Restore initial values */
  249. DBL_MEMCPY (y, y0, dim);
  250. return s;
  251. }
  252. }
  253. /* difference between 4th and 5th order */
  254. for (i = 0; i < dim; i++)
  255. {
  256. yerr[i] = h * (ec[1] * k1[i] + ec[3] * k3[i] + ec[4] * k4[i]
  257. + ec[5] * k5[i] + ec[6] * k6[i]);
  258. }
  259. return GSL_SUCCESS;
  260. }
  261. static int
  262. rkf45_reset (void *vstate, size_t dim)
  263. {
  264. rkf45_state_t *state = (rkf45_state_t *) vstate;
  265. DBL_ZERO_MEMSET (state->k1, dim);
  266. DBL_ZERO_MEMSET (state->k2, dim);
  267. DBL_ZERO_MEMSET (state->k3, dim);
  268. DBL_ZERO_MEMSET (state->k4, dim);
  269. DBL_ZERO_MEMSET (state->k5, dim);
  270. DBL_ZERO_MEMSET (state->k6, dim);
  271. DBL_ZERO_MEMSET (state->ytmp, dim);
  272. DBL_ZERO_MEMSET (state->y0, dim);
  273. return GSL_SUCCESS;
  274. }
  275. static unsigned int
  276. rkf45_order (void *vstate)
  277. {
  278. rkf45_state_t *state = (rkf45_state_t *) vstate;
  279. state = 0; /* prevent warnings about unused parameters */
  280. return 5;
  281. }
  282. static void
  283. rkf45_free (void *vstate)
  284. {
  285. rkf45_state_t *state = (rkf45_state_t *) vstate;
  286. free (state->ytmp);
  287. free (state->y0);
  288. free (state->k6);
  289. free (state->k5);
  290. free (state->k4);
  291. free (state->k3);
  292. free (state->k2);
  293. free (state->k1);
  294. free (state);
  295. }
  296. static const gsl_odeiv_step_type rkf45_type = { "rkf45", /* name */
  297. 1, /* can use dydt_in */
  298. 0, /* gives exact dydt_out */
  299. &rkf45_alloc,
  300. &rkf45_apply,
  301. &rkf45_reset,
  302. &rkf45_order,
  303. &rkf45_free
  304. };
  305. const gsl_odeiv_step_type *gsl_odeiv_step_rkf45 = &rkf45_type;