3 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
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.
27 #include <gsl/gsl_math.h>
28 #include <gsl/gsl_errno.h>
29 #include "odeiv_util.h"
30 #include <gsl/gsl_odeiv.h>
33 /* gear2 state object */
36 int primed; /* flag indicating that yim1 is ready */
37 double t_primed; /* system was primed for this value of t */
38 double last_h; /* last step size */
39 gsl_odeiv_step *primer; /* stepper to use for priming */
40 double *yim1; /* y_{i-1} */
41 double *k; /* work space */
42 double *y0; /* work space */
50 gear2_alloc (size_t dim)
52 gear2_state_t *state = (gear2_state_t *) malloc (sizeof (gear2_state_t));
56 GSL_ERROR_NULL ("failed to allocate space for gear2_state", GSL_ENOMEM);
59 state->yim1 = (double *) malloc (dim * sizeof (double));
64 GSL_ERROR_NULL ("failed to allocate space for yim1", GSL_ENOMEM);
67 state->k = (double *) malloc (dim * sizeof (double));
73 GSL_ERROR_NULL ("failed to allocate space for k", GSL_ENOMEM);
76 state->y0 = (double *) malloc (dim * sizeof (double));
83 GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
86 state->y0_orig = (double *) malloc (dim * sizeof (double));
88 if (state->y0_orig == 0)
94 GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
97 state->y_onestep = (double *) malloc (dim * sizeof (double));
99 if (state->y_onestep == 0)
101 free (state->y0_orig);
106 GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
110 state->primer = gsl_odeiv_step_alloc (gsl_odeiv_step_rk4imp, dim);
112 if (state->primer == 0)
114 free (state->y_onestep);
115 free (state->y0_orig);
120 GSL_ERROR_NULL ("failed to allocate space for primer", GSL_ENOMEM);
129 gear2_step (double *y, gear2_state_t * state,
130 const double h, const double t,
131 const size_t dim, const gsl_odeiv_system * sys)
133 /* Makes a Gear2 advance with step size h.
134 y0 is the initial values of variables y.
135 The implicit matrix equations to solve are:
136 k = y0 + h * f(t + h, k)
137 y = y0 + h * f(t + h, k)
140 const int iter_steps = 3;
143 double *y0 = state->y0;
144 double *yim1 = state->yim1;
145 double *k = state->k;
147 /* Iterative solution of k = y0 + h * f(t + h, k)
148 Note: This method does not check for convergence of the
152 for (nu = 0; nu < iter_steps; nu++)
154 int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, k);
156 if (s != GSL_SUCCESS)
161 for (i = 0; i < dim; i++)
163 y[i] = ((4.0 * y0[i] - yim1[i]) + 2.0 * h * k[i]) / 3.0;
171 gear2_apply (void *vstate,
177 const double dydt_in[],
178 double dydt_out[], const gsl_odeiv_system * sys)
180 gear2_state_t *state = (gear2_state_t *) vstate;
184 if (state->primed == 0 || t == state->t_primed || h != state->last_h)
186 /* Execute a single-step method to prime the process. Note that
187 * we do this if the step size changes, so frequent step size
188 * changes will cause the method to stutter.
190 * Note that we reuse this method if the time has not changed,
191 * which can occur when the adaptive driver is attempting to find
192 * an appropriate step-size on its first iteration */
195 DBL_MEMCPY (state->yim1, y, dim);
198 gsl_odeiv_step_apply (state->primer, t, h, y, yerr, dydt_in, dydt_out,
201 /* Make note of step size and indicate readiness for a Gear step. */
212 /* We have a previous y value in the buffer, and the step
213 * sizes match, so we go ahead with the Gear step.
216 double *const k = state->k;
217 double *const y0 = state->y0;
218 double *const y0_orig = state->y0_orig;
219 double *const yim1 = state->yim1;
220 double *y_onestep = state->y_onestep;
225 DBL_MEMCPY (y0, y, dim);
227 /* iterative solution */
229 if (dydt_out != NULL)
231 DBL_MEMCPY (k, dydt_out, dim);
234 /* First traverse h with one step (save to y_onestep) */
236 DBL_MEMCPY (y_onestep, y, dim);
238 s = gear2_step (y_onestep, state, h, t, dim, sys);
240 if (s != GSL_SUCCESS)
245 /* Then with two steps with half step length (save to y) */
247 s = gear2_step (y, state, h / 2.0, t, dim, sys);
249 if (s != GSL_SUCCESS)
251 /* Restore original y vector */
252 DBL_MEMCPY (y, y0_orig, dim);
256 DBL_MEMCPY (y0, y, dim);
258 s = gear2_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
260 if (s != GSL_SUCCESS)
262 /* Restore original y vector */
263 DBL_MEMCPY (y, y0_orig, dim);
269 if (dydt_out != NULL)
271 s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
273 if (s != GSL_SUCCESS)
275 /* Restore original y vector */
276 DBL_MEMCPY (y, y0_orig, dim);
281 /* Estimate error and update the state buffer. */
283 for (i = 0; i < dim; i++)
285 yerr[i] = 4.0 * (y[i] - y_onestep[i]);
289 /* Make note of step size. */
297 gear2_reset (void *vstate, size_t dim)
299 gear2_state_t *state = (gear2_state_t *) vstate;
301 DBL_ZERO_MEMSET (state->yim1, dim);
302 DBL_ZERO_MEMSET (state->k, dim);
303 DBL_ZERO_MEMSET (state->y0, dim);
311 gear2_order (void *vstate)
313 gear2_state_t *state = (gear2_state_t *) vstate;
314 state = 0; /* prevent warnings about unused parameters */
319 gear2_free (void *vstate)
321 gear2_state_t *state = (gear2_state_t *) vstate;
326 free (state->y0_orig);
327 free (state->y_onestep);
328 gsl_odeiv_step_free (state->primer);
333 static const gsl_odeiv_step_type gear2_type = { "gear2", /* name */
334 1, /* can use dydt_in */
335 0, /* gives exact dydt_out */
343 const gsl_odeiv_step_type *gsl_odeiv_step_gear2 = &gear2_type;