1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272
|
/* ode-initval/gear1.c
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
/* Gear 1. This is the implicit Euler a.k.a backward Euler method. */
/* Author: G. Jungman
*/
/* Error estimation by step doubling, see eg. Ascher, U.M., Petzold,
L.R., Computer methods for ordinary differential and
differential-algebraic equations, SIAM, Philadelphia, 1998.
The method is also described in eg. this reference.
*/
#include <config.h>
#include <stdlib.h>
#include <string.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_odeiv.h>
#include "odeiv_util.h"
typedef struct
{
double *k;
double *y0;
double *y0_orig;
double *y_onestep;
}
gear1_state_t;
static void *
gear1_alloc (size_t dim)
{
gear1_state_t *state = (gear1_state_t *) malloc (sizeof (gear1_state_t));
if (state == 0)
{
GSL_ERROR_NULL ("failed to allocate space for gear1_state", GSL_ENOMEM);
}
state->k = (double *) malloc (dim * sizeof (double));
if (state->k == 0)
{
free (state);
GSL_ERROR_NULL ("failed to allocate space for k", GSL_ENOMEM);
}
state->y0 = (double *) malloc (dim * sizeof (double));
if (state->y0 == 0)
{
free (state->k);
free (state);
GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
}
state->y0_orig = (double *) malloc (dim * sizeof (double));
if (state->y0_orig == 0)
{
free (state->y0);
free (state->k);
free (state);
GSL_ERROR_NULL ("failed to allocate space for y0_orig", GSL_ENOMEM);
}
state->y_onestep = (double *) malloc (dim * sizeof (double));
if (state->y_onestep == 0)
{
free (state->y0_orig);
free (state->y0);
free (state->k);
free (state);
GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
}
return state;
}
static int
gear1_step (double *y, gear1_state_t *state,
const double h, const double t,
const size_t dim, const gsl_odeiv_system *sys)
{
/* Makes an implicit Euler advance with step size h.
y0 is the initial values of variables y.
The implicit matrix equations to solve are:
k = y0 + h * f(t + h, k)
y = y0 + h * f(t + h, k)
*/
const int iter_steps = 3;
int nu;
size_t i;
double *y0 = state->y0;
double *k = state->k;
/* Iterative solution of k = y0 + h * f(t + h, k)
Note: This method does not check for convergence of the
iterative solution!
*/
for (nu = 0; nu < iter_steps; nu++)
{
int s = GSL_ODEIV_FN_EVAL(sys, t + h, y, k);
if (s != GSL_SUCCESS)
{
return s;
}
for (i=0; i<dim; i++)
{
y[i] = y0[i] + h * k[i];
}
}
return GSL_SUCCESS;
}
static int
gear1_apply(void * vstate,
size_t dim,
double t,
double h,
double y[],
double yerr[],
const double dydt_in[],
double dydt_out[],
const gsl_odeiv_system * sys)
{
gear1_state_t *state = (gear1_state_t *) vstate;
size_t i;
double *y0 = state->y0;
double *y0_orig = state->y0_orig;
double *y_onestep = state->y_onestep;
/* initialization */
DBL_MEMCPY(y0, y, dim);
/* Save initial values for possible failures */
DBL_MEMCPY (y0_orig, y, dim);
/* First traverse h with one step (save to y_onestep) */
DBL_MEMCPY (y_onestep, y, dim);
{
int s = gear1_step (y_onestep, state, h, t, dim, sys);
if (s != GSL_SUCCESS)
{
return s;
}
}
/* Then with two steps with half step length (save to y) */
{
int s = gear1_step (y, state, h / 2.0, t, dim, sys);
if (s != GSL_SUCCESS)
{
/* Restore original y vector */
DBL_MEMCPY (y, y0_orig, dim);
return s;
}
}
DBL_MEMCPY (y0, y, dim);
{
int s = gear1_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
if (s != GSL_SUCCESS)
{
/* Restore original y vector */
DBL_MEMCPY (y, y0_orig, dim);
return s;
}
}
/* Cleanup update */
if (dydt_out != NULL)
{
int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
if (s != GSL_SUCCESS)
{
/* Restore original y vector */
DBL_MEMCPY (y, y0_orig, dim);
return s;
}
}
/* Error estimation */
for (i = 0; i < dim; i++)
{
yerr[i] = 4.0 * (y[i] - y_onestep[i]);
}
return GSL_SUCCESS;
}
static int
gear1_reset (void *vstate, size_t dim)
{
gear1_state_t *state = (gear1_state_t *) vstate;
DBL_ZERO_MEMSET (state->y_onestep, dim);
DBL_ZERO_MEMSET (state->y0_orig, dim);
DBL_ZERO_MEMSET (state->y0, dim);
DBL_ZERO_MEMSET (state->k, dim);
return GSL_SUCCESS;
}
static unsigned int
gear1_order (void *vstate)
{
gear1_state_t *state = (gear1_state_t *) vstate;
state = 0; /* prevent warnings about unused parameters */
return 1;
}
static void
gear1_free (void *vstate)
{
gear1_state_t *state = (gear1_state_t *) vstate;
free (state->y_onestep);
free (state->y0_orig);
free (state->y0);
free (state->k);
free (state);
}
static const gsl_odeiv_step_type gear1_type = { "gear1", /* name */
0, /* can use dydt_in? */
1, /* gives exact dydt_out? */
&gear1_alloc,
&gear1_apply,
&gear1_reset,
&gear1_order,
&gear1_free
};
const gsl_odeiv_step_type *gsl_odeiv_step_gear1 = &gear1_type;
|