Blame ode-initval/gear1.c

Packit 67cb25
/* ode-initval/gear1.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
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
/* Gear 1. This is the implicit Euler a.k.a backward Euler method. */
Packit 67cb25
Packit 67cb25
/* Author:  G. Jungman
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
/* Error estimation by step doubling, see eg. Ascher, U.M., Petzold,
Packit 67cb25
   L.R., Computer methods for ordinary differential and
Packit 67cb25
   differential-algebraic equations, SIAM, Philadelphia, 1998.
Packit 67cb25
   The method is also described in eg. this reference.
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_odeiv.h>
Packit 67cb25
Packit 67cb25
#include "odeiv_util.h"
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  double *k;
Packit 67cb25
  double *y0;
Packit 67cb25
  double *y0_orig;
Packit 67cb25
  double *y_onestep;
Packit 67cb25
}
Packit 67cb25
gear1_state_t;
Packit 67cb25
Packit 67cb25
static void *
Packit 67cb25
gear1_alloc (size_t dim)
Packit 67cb25
{
Packit 67cb25
  gear1_state_t *state = (gear1_state_t *) malloc (sizeof (gear1_state_t));
Packit 67cb25
Packit 67cb25
  if (state == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for gear1_state", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->y0 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->y0 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for y0", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->y0_orig = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->y0_orig == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->y0);
Packit 67cb25
      free (state->k);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for y0_orig", 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
      free (state->y0_orig);
Packit 67cb25
      free (state->y0);
Packit 67cb25
      free (state->k);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for y_onestep", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return state;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
gear1_step (double *y, gear1_state_t *state, 
Packit 67cb25
	    const double h, const double t, 
Packit 67cb25
	    const size_t dim, const gsl_odeiv_system *sys)
Packit 67cb25
{
Packit 67cb25
  /* Makes an implicit Euler advance with step size h.
Packit 67cb25
     y0 is the initial values of variables y. 
Packit 67cb25
Packit 67cb25
     The implicit matrix equations to solve are:
Packit 67cb25
Packit 67cb25
     k = y0 + h * f(t + h, k)
Packit 67cb25
Packit 67cb25
     y = y0 + h * f(t + h, k)
Packit 67cb25
  */
Packit 67cb25
Packit 67cb25
  const int iter_steps = 3;
Packit 67cb25
  int nu;
Packit 67cb25
  size_t i;
Packit 67cb25
  double *y0 = state->y0;
Packit 67cb25
  double *k = state->k;
Packit 67cb25
Packit 67cb25
  /* Iterative solution of k = y0 + h * f(t + h, k)
Packit 67cb25
Packit 67cb25
     Note: This method does not check for convergence of the
Packit 67cb25
     iterative solution! 
Packit 67cb25
  */
Packit 67cb25
Packit 67cb25
  for (nu = 0; nu < iter_steps; nu++) 
Packit 67cb25
    {
Packit 67cb25
      int s = GSL_ODEIV_FN_EVAL(sys, t + h, y, k);
Packit 67cb25
Packit 67cb25
      if (s != GSL_SUCCESS)
Packit 67cb25
	{
Packit 67cb25
	  return s;
Packit 67cb25
	}     
Packit 67cb25
Packit 67cb25
      for (i=0; i
Packit 67cb25
	{
Packit 67cb25
	  y[i] = y0[i] + h * k[i];
Packit 67cb25
	}
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
gear1_apply(void * vstate,
Packit 67cb25
            size_t dim,
Packit 67cb25
            double t,
Packit 67cb25
            double h,
Packit 67cb25
            double y[],
Packit 67cb25
            double yerr[],
Packit 67cb25
            const double dydt_in[],
Packit 67cb25
            double dydt_out[],
Packit 67cb25
            const gsl_odeiv_system * sys)
Packit 67cb25
{
Packit 67cb25
  gear1_state_t *state = (gear1_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  double *y0 = state->y0;
Packit 67cb25
  double *y0_orig = state->y0_orig;
Packit 67cb25
  double *y_onestep = state->y_onestep;
Packit 67cb25
Packit 67cb25
  /* initialization */
Packit 67cb25
  DBL_MEMCPY(y0, y, dim);
Packit 67cb25
Packit 67cb25
  /* Save initial values for possible failures */
Packit 67cb25
  DBL_MEMCPY (y0_orig, y, dim);
Packit 67cb25
Packit 67cb25
  /* First traverse h with one step (save to y_onestep) */
Packit 67cb25
  DBL_MEMCPY (y_onestep, y, dim);
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int s = gear1_step (y_onestep, state, h, t, dim, sys);
Packit 67cb25
Packit 67cb25
    if (s != GSL_SUCCESS) 
Packit 67cb25
      {
Packit 67cb25
        return s;
Packit 67cb25
      }
Packit 67cb25
  }    
Packit 67cb25
Packit 67cb25
 /* Then with two steps with half step length (save to y) */ 
Packit 67cb25
  {
Packit 67cb25
    int s = gear1_step (y, state, h / 2.0, t, dim, sys);
Packit 67cb25
Packit 67cb25
    if (s != GSL_SUCCESS) 
Packit 67cb25
      {
Packit 67cb25
        /* Restore original y vector */
Packit 67cb25
        DBL_MEMCPY (y, y0_orig, dim);
Packit 67cb25
        return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  DBL_MEMCPY (y0, y, dim);
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int s = gear1_step (y, state, h / 2.0, t + h / 2.0, dim, sys);
Packit 67cb25
Packit 67cb25
    if (s != GSL_SUCCESS) 
Packit 67cb25
      {
Packit 67cb25
        /* Restore original y vector */
Packit 67cb25
        DBL_MEMCPY (y, y0_orig, dim);
Packit 67cb25
        return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
  
Packit 67cb25
  /* Cleanup update */
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 y vector */
Packit 67cb25
          DBL_MEMCPY (y, y0_orig, dim);
Packit 67cb25
          return s;
Packit 67cb25
        } 
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  /* Error estimation */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++) 
Packit 67cb25
    {
Packit 67cb25
      yerr[i] = 4.0 * (y[i] - y_onestep[i]);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
gear1_reset (void *vstate, size_t dim)
Packit 67cb25
{
Packit 67cb25
  gear1_state_t *state = (gear1_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  DBL_ZERO_MEMSET (state->y_onestep, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->y0_orig, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->y0, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k, dim);
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static unsigned int
Packit 67cb25
gear1_order (void *vstate)
Packit 67cb25
{
Packit 67cb25
  gear1_state_t *state = (gear1_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
gear1_free (void *vstate)
Packit 67cb25
{
Packit 67cb25
  gear1_state_t *state = (gear1_state_t *) vstate;
Packit 67cb25
  free (state->y_onestep);
Packit 67cb25
  free (state->y0_orig);
Packit 67cb25
  free (state->y0);
Packit 67cb25
  free (state->k);
Packit 67cb25
  free (state);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static const gsl_odeiv_step_type gear1_type = { "gear1",        /* name */
Packit 67cb25
  0,                            /* can use dydt_in? */
Packit 67cb25
  1,                            /* gives exact dydt_out? */
Packit 67cb25
  &gear1_alloc,
Packit 67cb25
  &gear1_apply,
Packit 67cb25
  &gear1_reset,
Packit 67cb25
  &gear1_order,
Packit 67cb25
  &gear1_free
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
const gsl_odeiv_step_type *gsl_odeiv_step_gear1 = &gear1_type;