3 * Copyright (C) 2001, 2004, 2007 Brian Gough
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.
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.
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.
20 /* Runge-Kutta-Fehlberg 4(5)*/
22 /* Reference eg. Hairer, E., Norsett S.P., Wanner, G. Solving ordinary
23 differential equations I, Nonstiff Problems, 2nd revised edition,
30 #include <gsl/gsl_errno.h>
31 #include <gsl/gsl_odeiv.h>
33 #include "odeiv_util.h"
35 /* Runge-Kutta-Fehlberg coefficients. Zero elements left out */
37 static const double ah[] = { 1.0/4.0, 3.0/8.0, 12.0/13.0, 1.0, 1.0/2.0 };
38 static const double b3[] = { 3.0/32.0, 9.0/32.0 };
39 static const double b4[] = { 1932.0/2197.0, -7200.0/2197.0, 7296.0/2197.0};
40 static const double b5[] = { 8341.0/4104.0, -32832.0/4104.0, 29440.0/4104.0, -845.0/4104.0};
41 static const double b6[] = { -6080.0/20520.0, 41040.0/20520.0, -28352.0/20520.0, 9295.0/20520.0, -5643.0/20520.0};
43 static const double c1 = 902880.0/7618050.0;
44 static const double c3 = 3953664.0/7618050.0;
45 static const double c4 = 3855735.0/7618050.0;
46 static const double c5 = -1371249.0/7618050.0;
47 static const double c6 = 277020.0/7618050.0;
49 /* These are the differences of fifth and fourth order coefficients
50 for error estimation */
52 static const double ec[] = { 0.0,
75 rkf45_alloc (size_t dim)
77 rkf45_state_t *state = (rkf45_state_t *) malloc (sizeof (rkf45_state_t));
81 GSL_ERROR_NULL ("failed to allocate space for rkf45_state", GSL_ENOMEM);
84 state->k1 = (double *) malloc (dim * sizeof (double));
89 GSL_ERROR_NULL ("failed to allocate space for k1", GSL_ENOMEM);
92 state->k2 = (double *) malloc (dim * sizeof (double));
98 GSL_ERROR_NULL ("failed to allocate space for k2", GSL_ENOMEM);
101 state->k3 = (double *) malloc (dim * sizeof (double));
108 GSL_ERROR_NULL ("failed to allocate space for k3", GSL_ENOMEM);
111 state->k4 = (double *) malloc (dim * sizeof (double));
119 GSL_ERROR_NULL ("failed to allocate space for k4", GSL_ENOMEM);
122 state->k5 = (double *) malloc (dim * sizeof (double));
131 GSL_ERROR_NULL ("failed to allocate space for k5", GSL_ENOMEM);
134 state->k6 = (double *) malloc (dim * sizeof (double));
144 GSL_ERROR_NULL ("failed to allocate space for k6", GSL_ENOMEM);
147 state->y0 = (double *) malloc (dim * sizeof (double));
158 GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
161 state->ytmp = (double *) malloc (dim * sizeof (double));
163 if (state->ytmp == 0)
173 GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
181 rkf45_apply (void *vstate,
187 const double dydt_in[],
188 double dydt_out[], const gsl_odeiv_system * sys)
190 rkf45_state_t *state = (rkf45_state_t *) vstate;
194 double *const k1 = state->k1;
195 double *const k2 = state->k2;
196 double *const k3 = state->k3;
197 double *const k4 = state->k4;
198 double *const k5 = state->k5;
199 double *const k6 = state->k6;
200 double *const ytmp = state->ytmp;
201 double *const y0 = state->y0;
203 DBL_MEMCPY (y0, y, dim);
208 DBL_MEMCPY (k1, dydt_in, dim);
212 int s = GSL_ODEIV_FN_EVAL (sys, t, y, k1);
214 if (s != GSL_SUCCESS)
220 for (i = 0; i < dim; i++)
221 ytmp[i] = y[i] + ah[0] * h * k1[i];
225 int s = GSL_ODEIV_FN_EVAL (sys, t + ah[0] * h, ytmp, k2);
227 if (s != GSL_SUCCESS)
233 for (i = 0; i < dim; i++)
234 ytmp[i] = y[i] + h * (b3[0] * k1[i] + b3[1] * k2[i]);
238 int s = GSL_ODEIV_FN_EVAL (sys, t + ah[1] * h, ytmp, k3);
240 if (s != GSL_SUCCESS)
246 for (i = 0; i < dim; i++)
247 ytmp[i] = y[i] + h * (b4[0] * k1[i] + b4[1] * k2[i] + b4[2] * k3[i]);
251 int s = GSL_ODEIV_FN_EVAL (sys, t + ah[2] * h, ytmp, k4);
253 if (s != GSL_SUCCESS)
259 for (i = 0; i < dim; i++)
261 y[i] + h * (b5[0] * k1[i] + b5[1] * k2[i] + b5[2] * k3[i] +
266 int s = GSL_ODEIV_FN_EVAL (sys, t + ah[3] * h, ytmp, k5);
268 if (s != GSL_SUCCESS)
274 for (i = 0; i < dim; i++)
276 y[i] + h * (b6[0] * k1[i] + b6[1] * k2[i] + b6[2] * k3[i] +
277 b6[3] * k4[i] + b6[4] * k5[i]);
279 /* k6 step and final sum */
281 int s = GSL_ODEIV_FN_EVAL (sys, t + ah[4] * h, ytmp, k6);
283 if (s != GSL_SUCCESS)
289 for (i = 0; i < dim; i++)
291 const double d_i = c1 * k1[i] + c3 * k3[i] + c4 * k4[i] + c5 * k5[i] + c6 * k6[i];
295 /* Derivatives at output */
297 if (dydt_out != NULL)
299 int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
301 if (s != GSL_SUCCESS)
303 /* Restore initial values */
304 DBL_MEMCPY (y, y0, dim);
310 /* difference between 4th and 5th order */
311 for (i = 0; i < dim; i++)
313 yerr[i] = h * (ec[1] * k1[i] + ec[3] * k3[i] + ec[4] * k4[i]
314 + ec[5] * k5[i] + ec[6] * k6[i]);
322 rkf45_reset (void *vstate, size_t dim)
324 rkf45_state_t *state = (rkf45_state_t *) vstate;
326 DBL_ZERO_MEMSET (state->k1, dim);
327 DBL_ZERO_MEMSET (state->k2, dim);
328 DBL_ZERO_MEMSET (state->k3, dim);
329 DBL_ZERO_MEMSET (state->k4, dim);
330 DBL_ZERO_MEMSET (state->k5, dim);
331 DBL_ZERO_MEMSET (state->k6, dim);
332 DBL_ZERO_MEMSET (state->ytmp, dim);
333 DBL_ZERO_MEMSET (state->y0, dim);
339 rkf45_order (void *vstate)
341 rkf45_state_t *state = (rkf45_state_t *) vstate;
342 state = 0; /* prevent warnings about unused parameters */
347 rkf45_free (void *vstate)
349 rkf45_state_t *state = (rkf45_state_t *) vstate;
362 static const gsl_odeiv_step_type rkf45_type = { "rkf45", /* name */
363 1, /* can use dydt_in */
364 0, /* gives exact dydt_out */
372 const gsl_odeiv_step_type *gsl_odeiv_step_rkf45 = &rkf45_type;