|
Packit |
67cb25 |
/* ode-initval2/rk1imp.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2009, 2010 Tuomo Keskitalo
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is free software; you can redistribute it and/or modify
|
|
Packit |
67cb25 |
* it under the terms of the GNU General Public License as published by
|
|
Packit |
67cb25 |
* the Free Software Foundation; either version 3 of the License, or (at
|
|
Packit |
67cb25 |
* your option) any later version.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is distributed in the hope that it will be useful, but
|
|
Packit |
67cb25 |
* WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
Packit |
67cb25 |
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
Packit |
67cb25 |
* General Public License for more details.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* You should have received a copy of the GNU General Public License
|
|
Packit |
67cb25 |
* along with this program; if not, write to the Free Software
|
|
Packit |
67cb25 |
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Implicit Euler a.k.a backward Euler method. */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Reference: Ascher, U.M., Petzold, L.R., Computer methods for
|
|
Packit |
67cb25 |
ordinary differential and differential-algebraic equations, SIAM,
|
|
Packit |
67cb25 |
Philadelphia, 1998.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include <config.h>
|
|
Packit |
67cb25 |
#include <stdlib.h>
|
|
Packit |
67cb25 |
#include <string.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_math.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_errno.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_odeiv2.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include "odeiv_util.h"
|
|
Packit |
67cb25 |
#include "rksubs.c"
|
|
Packit |
67cb25 |
#include "modnewton1.c"
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Stage of method */
|
|
Packit |
67cb25 |
#define RK1IMP_STAGE 1
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_matrix *A; /* Runge-Kutta coefficients */
|
|
Packit |
67cb25 |
double *y_onestep; /* Result with one step */
|
|
Packit |
67cb25 |
double *y_twostep; /* Result with two half steps */
|
|
Packit |
67cb25 |
double *ytmp; /* Temporary work space */
|
|
Packit |
67cb25 |
double *y_save; /* Backup space */
|
|
Packit |
67cb25 |
double *YZ; /* Runge-Kutta points */
|
|
Packit |
67cb25 |
double *fYZ; /* Derivatives at YZ */
|
|
Packit |
67cb25 |
gsl_matrix *dfdy; /* Jacobian matrix */
|
|
Packit |
67cb25 |
double *dfdt; /* time derivative of f */
|
|
Packit |
67cb25 |
modnewton1_state_t *esol; /* nonlinear equation solver */
|
|
Packit |
67cb25 |
double *errlev; /* desired error level of y */
|
|
Packit |
67cb25 |
const gsl_odeiv2_driver *driver; /* pointer to driver object */
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
rk1imp_state_t;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void *
|
|
Packit |
67cb25 |
rk1imp_alloc (size_t dim)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) malloc (sizeof (rk1imp_state_t));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for rk1imp_state",
|
|
Packit |
67cb25 |
GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->A = gsl_matrix_alloc (RK1IMP_STAGE, RK1IMP_STAGE);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->A == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for A", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->y_onestep = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->y_onestep == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->y_twostep = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->y_twostep == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->ytmp = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->ytmp == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->y_save = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->y_save == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for y_save", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->YZ = (double *) malloc (dim * RK1IMP_STAGE * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->YZ == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for YZ", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->fYZ = (double *) malloc (dim * RK1IMP_STAGE * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->fYZ == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for fYZ", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->dfdt = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->dfdt == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->fYZ);
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for dfdt", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->dfdy = gsl_matrix_alloc (dim, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->dfdy == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->fYZ);
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for dfdy", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->esol = modnewton1_alloc (dim, RK1IMP_STAGE);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->esol == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_matrix_free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->fYZ);
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for esol", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->errlev = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->errlev == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
modnewton1_free (state->esol);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->fYZ);
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for errlev", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->driver = NULL;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return state;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk1imp_apply (void *vstate, size_t dim, double t, double h,
|
|
Packit |
67cb25 |
double y[], double yerr[],
|
|
Packit |
67cb25 |
const double dydt_in[], double dydt_out[],
|
|
Packit |
67cb25 |
const gsl_odeiv2_system * sys)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* Makes an implicit Euler step with size h and estimates the local
|
|
Packit |
67cb25 |
error of the step by step doubling.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double *const y_onestep = state->y_onestep;
|
|
Packit |
67cb25 |
double *const y_twostep = state->y_twostep;
|
|
Packit |
67cb25 |
double *const ytmp = state->ytmp;
|
|
Packit |
67cb25 |
double *const y_save = state->y_save;
|
|
Packit |
67cb25 |
double *const YZ = state->YZ;
|
|
Packit |
67cb25 |
double *const fYZ = state->fYZ;
|
|
Packit |
67cb25 |
gsl_matrix *const dfdy = state->dfdy;
|
|
Packit |
67cb25 |
double *const dfdt = state->dfdt;
|
|
Packit |
67cb25 |
double *const errlev = state->errlev;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const modnewton1_state_t *esol = state->esol;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Runge-Kutta coefficients */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix *A = state->A;
|
|
Packit |
67cb25 |
const double b[] = { 1.0 };
|
|
Packit |
67cb25 |
const double c[] = { 1.0 };
|
|
Packit |
67cb25 |
gsl_matrix_set (A, 0, 0, 1.0);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (esol == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("no non-linear equation solver speficied", GSL_EINVAL);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Get desired error levels via gsl_odeiv2_control object through driver
|
|
Packit |
67cb25 |
object, which is a requirement for this stepper.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->driver == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return GSL_EFAULT;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < dim; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if (dydt_in != NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_odeiv2_control_errlevel (state->driver->c, y[i],
|
|
Packit |
67cb25 |
dydt_in[i], h, i, &errlev[i]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_odeiv2_control_errlevel (state->driver->c, y[i],
|
|
Packit |
67cb25 |
0.0, h, i, &errlev[i]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Evaluate Jacobian for modnewton1 */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = GSL_ODEIV_JA_EVAL (sys, t, y, dfdy->data, dfdt);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Calculate a single step with size h */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = modnewton1_init ((void *) esol, A, h, dfdy, sys);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = modnewton1_solve ((void *) esol, A, c, t, h, y,
|
|
Packit |
67cb25 |
sys, YZ, errlev);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = GSL_ODEIV_FN_EVAL (sys, t + c[0] * h, YZ, fYZ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = rksubs (y_onestep, h, y, fYZ, b, RK1IMP_STAGE, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Error estimation by step doubling */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = modnewton1_init ((void *) esol, A, h / 2.0, dfdy, sys);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* 1st half step */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = modnewton1_solve ((void *) esol, A, c, t, h / 2.0, y,
|
|
Packit |
67cb25 |
sys, YZ, errlev);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = GSL_ODEIV_FN_EVAL (sys, t + c[0] * h / 2.0, YZ, fYZ);
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = rksubs (ytmp, h / 2.0, y, fYZ, b, RK1IMP_STAGE, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Save original y values in case of error */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
DBL_MEMCPY (y_save, y, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* 2nd half step */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = modnewton1_solve ((void *) esol, A, c, t + h / 2.0, h / 2.0,
|
|
Packit |
67cb25 |
ytmp, sys, YZ, errlev);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = GSL_ODEIV_FN_EVAL (sys, t + h / 2.0 + c[0] * h / 2.0, YZ, fYZ);
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* Note: rk1imp returns y using the results from two half steps
|
|
Packit |
67cb25 |
instead of the single step since the results are freely
|
|
Packit |
67cb25 |
available and more precise.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int s = rksubs (y_twostep, h / 2.0, ytmp, fYZ, b, RK1IMP_STAGE, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
DBL_MEMCPY (y, y_save, dim);
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
DBL_MEMCPY (y, y_twostep, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Error estimation */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
for (i = 0; i < dim; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
yerr[i] = ODEIV_ERR_SAFETY * 0.5 * fabs (y_twostep[i] - y_onestep[i]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Derivatives at output */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (dydt_out != NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = GSL_ODEIV_FN_EVAL (sys, t + h, y, dydt_out);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* Restore original values */
|
|
Packit |
67cb25 |
DBL_MEMCPY (y, y_save, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk1imp_set_driver (void *vstate, const gsl_odeiv2_driver * d)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->driver = d;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk1imp_reset (void *vstate, size_t dim)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y_onestep, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y_twostep, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->ytmp, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y_save, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->YZ, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->fYZ, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static unsigned int
|
|
Packit |
67cb25 |
rk1imp_order (void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) vstate;
|
|
Packit |
67cb25 |
state = 0; /* prevent warnings about unused parameters */
|
|
Packit |
67cb25 |
return 1;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
rk1imp_free (void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk1imp_state_t *state = (rk1imp_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free (state->errlev);
|
|
Packit |
67cb25 |
modnewton1_free (state->esol);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->fYZ);
|
|
Packit |
67cb25 |
free (state->YZ);
|
|
Packit |
67cb25 |
free (state->y_save);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->y_twostep);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_matrix_free (state->A);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static const gsl_odeiv2_step_type rk1imp_type = {
|
|
Packit |
67cb25 |
"rk1imp", /* name */
|
|
Packit |
67cb25 |
1, /* can use dydt_in? */
|
|
Packit |
67cb25 |
1, /* gives exact dydt_out? */
|
|
Packit |
67cb25 |
&rk1imp_alloc,
|
|
Packit |
67cb25 |
&rk1imp_apply,
|
|
Packit |
67cb25 |
&rk1imp_set_driver,
|
|
Packit |
67cb25 |
&rk1imp_reset,
|
|
Packit |
67cb25 |
&rk1imp_order,
|
|
Packit |
67cb25 |
&rk1imp_free
|
|
Packit |
67cb25 |
};
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const gsl_odeiv2_step_type *gsl_odeiv2_step_rk1imp = &rk1imp_type;
|