Blame ode-initval/rk4imp.c

Packit 67cb25
/* ode-initval/rk4imp.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
/* Runge-Kutta 4, Gaussian implicit */
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
   Method coefficients can also be found in it.
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 *k1nu;
Packit 67cb25
  double *k2nu;
Packit 67cb25
  double *ytmp1;
Packit 67cb25
  double *ytmp2;
Packit 67cb25
  double *y0;
Packit 67cb25
  double *y0_orig;
Packit 67cb25
  double *y_onestep;
Packit 67cb25
}
Packit 67cb25
rk4imp_state_t;
Packit 67cb25
Packit 67cb25
static void *
Packit 67cb25
rk4imp_alloc (size_t dim)
Packit 67cb25
{
Packit 67cb25
  rk4imp_state_t *state = (rk4imp_state_t *) malloc (sizeof (rk4imp_state_t));
Packit 67cb25
Packit 67cb25
  if (state == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for rk4imp_state",
Packit 67cb25
                      GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k1nu = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k1nu == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k1nu", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k2nu = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k2nu == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k1nu);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k2nu", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->ytmp1 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->ytmp1 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k2nu);
Packit 67cb25
      free (state->k1nu);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for ytmp1", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->ytmp2 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->ytmp2 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->ytmp1);
Packit 67cb25
      free (state->k2nu);
Packit 67cb25
      free (state->k1nu);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for ytmp2", 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->ytmp2);
Packit 67cb25
      free (state->ytmp1);
Packit 67cb25
      free (state->k2nu);
Packit 67cb25
      free (state->k1nu);
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->ytmp2);
Packit 67cb25
      free (state->ytmp1);
Packit 67cb25
      free (state->k2nu);
Packit 67cb25
      free (state->k1nu);
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->ytmp2);
Packit 67cb25
      free (state->ytmp1);
Packit 67cb25
      free (state->k2nu);
Packit 67cb25
      free (state->k1nu);
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
rk4imp_step (double *y, rk4imp_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 a Runge-Kutta 4th order implicit advance with step size h.
Packit 67cb25
     y0 is initial values of variables y. 
Packit 67cb25
Packit 67cb25
     The implicit matrix equations to solve are:
Packit 67cb25
Packit 67cb25
     Y1 = y0 + h * a11 * f(t + h * c1, Y1) + h * a12 * f(t + h * c2, Y2) 
Packit 67cb25
     Y2 = y0 + h * a21 * f(t + h * c1, Y1) + h * a22 * f(t + h * c2, Y2) 
Packit 67cb25
Packit 67cb25
     y = y0 + h * b1 * f(t + h * c1, Y1) + h * b2 * f(t + h * c2, Y2)
Packit 67cb25
Packit 67cb25
     with constant coefficients a, b and c. For this method
Packit 67cb25
     they are: b=[0.5 0.5] c=[(3-sqrt(3))/6 (3+sqrt(3))/6]
Packit 67cb25
     a11=1/4, a12=(3-2*sqrt(3))/12, a21=(3+2*sqrt(3))/12 and a22=1/4
Packit 67cb25
  */
Packit 67cb25
Packit 67cb25
  const double ir3 = 1.0 / M_SQRT3;
Packit 67cb25
  const int iter_steps = 3;
Packit 67cb25
  int nu;
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  double *const k1nu = state->k1nu;
Packit 67cb25
  double *const k2nu = state->k2nu;
Packit 67cb25
  double *const ytmp1 = state->ytmp1;
Packit 67cb25
  double *const ytmp2 = state->ytmp2;
Packit 67cb25
Packit 67cb25
  /* iterative solution of Y1 and Y2.
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
      for (i = 0; i < dim; i++)
Packit 67cb25
        {
Packit 67cb25
          ytmp1[i] =
Packit 67cb25
            y[i] + h * (0.25 * k1nu[i] + 0.5 * (0.5 - ir3) * k2nu[i]);
Packit 67cb25
          ytmp2[i] =
Packit 67cb25
            y[i] + h * (0.25 * k2nu[i] + 0.5 * (0.5 + ir3) * k1nu[i]);
Packit 67cb25
        }
Packit 67cb25
      {
Packit 67cb25
        int s =
Packit 67cb25
	  GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h * (1.0 - ir3), ytmp1, k1nu);
Packit 67cb25
	
Packit 67cb25
	if (s != GSL_SUCCESS)
Packit 67cb25
	  {
Packit 67cb25
	    return s;
Packit 67cb25
	  }    
Packit 67cb25
      }
Packit 67cb25
      {
Packit 67cb25
        int s =
Packit 67cb25
	  GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h * (1.0 + ir3), ytmp2, k2nu);
Packit 67cb25
	
Packit 67cb25
	if (s != GSL_SUCCESS)
Packit 67cb25
	  {
Packit 67cb25
	    return s;
Packit 67cb25
	  }    
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  /* assignment */
Packit 67cb25
  
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      const double d_i = 0.5 * (k1nu[i] + k2nu[i]);
Packit 67cb25
      y[i] += h * d_i;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rk4imp_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
  rk4imp_state_t *state = (rk4imp_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
  double *k1nu = state->k1nu;
Packit 67cb25
  double *k2nu = state->k2nu;
Packit 67cb25
Packit 67cb25
  /* Initialization step */
Packit 67cb25
  DBL_MEMCPY (y0, y, dim);
Packit 67cb25
Packit 67cb25
  /* Save initial values in case of failure */
Packit 67cb25
  DBL_MEMCPY (y0_orig, y, dim);
Packit 67cb25
Packit 67cb25
  if (dydt_in != 0)
Packit 67cb25
    {
Packit 67cb25
      DBL_MEMCPY (k1nu, dydt_in, dim);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int s = GSL_ODEIV_FN_EVAL (sys, t, y, k1nu);
Packit 67cb25
Packit 67cb25
      if (s != GSL_SUCCESS)
Packit 67cb25
        {
Packit 67cb25
          return s;
Packit 67cb25
        }
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  DBL_MEMCPY (k2nu, k1nu, dim);
Packit 67cb25
Packit 67cb25
  /* First traverse h with one step (save to y_onestep) */
Packit 67cb25
Packit 67cb25
  DBL_MEMCPY (y_onestep, y, dim);
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int s = rk4imp_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
  {
Packit 67cb25
    int s = rk4imp_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 = GSL_ODEIV_FN_EVAL (sys, t + h/2.0, y, k1nu);
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 (k2nu, k1nu, dim);
Packit 67cb25
  
Packit 67cb25
  {
Packit 67cb25
    int s = rk4imp_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
  /* 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
	/* 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
  /* Denominator in step doubling error equation 
Packit 67cb25
   *  yerr = 0.5 * | y(onestep) - y(twosteps) | / (2^order - 1)  
Packit 67cb25
   */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++) 
Packit 67cb25
    {
Packit 67cb25
      yerr[i] = 8.0 * 0.5 * (y[i] - y_onestep[i]) / 15.0;
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rk4imp_reset (void *vstate, size_t dim)
Packit 67cb25
{
Packit 67cb25
  rk4imp_state_t *state = (rk4imp_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->k1nu, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k2nu, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->ytmp1, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->ytmp2, dim);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static unsigned int
Packit 67cb25
rk4imp_order (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rk4imp_state_t *state = (rk4imp_state_t *) vstate;
Packit 67cb25
  state = 0; /* prevent warnings about unused parameters */
Packit 67cb25
  return 4;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
rk4imp_free (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rk4imp_state_t *state = (rk4imp_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  free (state->y_onestep);
Packit 67cb25
  free (state->y0_orig);
Packit 67cb25
  free (state->y0);
Packit 67cb25
  free (state->k1nu);
Packit 67cb25
  free (state->k2nu);
Packit 67cb25
  free (state->ytmp1);
Packit 67cb25
  free (state->ytmp2);
Packit 67cb25
  free (state);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static const gsl_odeiv_step_type rk4imp_type = { "rk4imp",      /* name */
Packit 67cb25
  1,                             /* can use dydt_in? */
Packit 67cb25
  1,                             /* gives exact dydt_out? */
Packit 67cb25
  &rk4imp_alloc,
Packit 67cb25
  &rk4imp_apply,
Packit 67cb25
  &rk4imp_reset,
Packit 67cb25
  &rk4imp_order,
Packit 67cb25
  &rk4imp_free
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
const gsl_odeiv_step_type *gsl_odeiv_step_rk4imp = &rk4imp_type;