Blame ode-initval/rkck.c

Packit 67cb25
/* ode-initval/rkck.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(5), Cash-Karp */
Packit 67cb25
Packit 67cb25
/* Reference: Cash, J.R., Karp, A.H., ACM Transactions of Mathematical
Packit 67cb25
   Software, vol. 16 (1990) 201-222. */
Packit 67cb25
Packit 67cb25
/* Author:  G. Jungman
Packit 67cb25
 */
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <stdlib.h>
Packit 67cb25
#include <string.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
/* Cash-Karp constants */
Packit 67cb25
static const double ah[] = { 1.0 / 5.0, 0.3, 3.0 / 5.0, 1.0, 7.0 / 8.0 };
Packit 67cb25
static const double b21 = 1.0 / 5.0;
Packit 67cb25
static const double b3[] = { 3.0 / 40.0, 9.0 / 40.0 };
Packit 67cb25
static const double b4[] = { 0.3, -0.9, 1.2 };
Packit 67cb25
static const double b5[] = { -11.0 / 54.0, 2.5, -70.0 / 27.0, 35.0 / 27.0 };
Packit 67cb25
static const double b6[] =
Packit 67cb25
  { 1631.0 / 55296.0, 175.0 / 512.0, 575.0 / 13824.0, 44275.0 / 110592.0,
Packit 67cb25
    253.0 / 4096.0 };
Packit 67cb25
static const double c1 = 37.0 / 378.0;
Packit 67cb25
static const double c3 = 250.0 / 621.0;
Packit 67cb25
static const double c4 = 125.0 / 594.0;
Packit 67cb25
static const double c6 = 512.0 / 1771.0;
Packit 67cb25
Packit 67cb25
/* These are the differences of fifth and fourth order coefficients
Packit 67cb25
   for error estimation */
Packit 67cb25
Packit 67cb25
static const double ec[] = { 0.0,
Packit 67cb25
  37.0 / 378.0 - 2825.0 / 27648.0,
Packit 67cb25
  0.0,
Packit 67cb25
  250.0 / 621.0 - 18575.0 / 48384.0,
Packit 67cb25
  125.0 / 594.0 - 13525.0 / 55296.0,
Packit 67cb25
  -277.0 / 14336.0,
Packit 67cb25
  512.0 / 1771.0 - 0.25
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  double *k1;
Packit 67cb25
  double *k2;
Packit 67cb25
  double *k3;
Packit 67cb25
  double *k4;
Packit 67cb25
  double *k5;
Packit 67cb25
  double *k6;
Packit 67cb25
  double *y0;
Packit 67cb25
  double *ytmp;
Packit 67cb25
}
Packit 67cb25
rkck_state_t;
Packit 67cb25
Packit 67cb25
static void *
Packit 67cb25
rkck_alloc (size_t dim)
Packit 67cb25
{
Packit 67cb25
  rkck_state_t *state = (rkck_state_t *) malloc (sizeof (rkck_state_t));
Packit 67cb25
Packit 67cb25
  if (state == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for rkck_state", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k1 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k1 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k1", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k2 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k2 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k2", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k3 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k3 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k3", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k4 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k4 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k3);
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k4", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k5 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k5 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k4);
Packit 67cb25
      free (state->k3);
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k5", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->k6 = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->k6 == 0)
Packit 67cb25
    {
Packit 67cb25
      free (state->k5);
Packit 67cb25
      free (state->k4);
Packit 67cb25
      free (state->k3);
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for k6", 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->k6);
Packit 67cb25
      free (state->k5);
Packit 67cb25
      free (state->k4);
Packit 67cb25
      free (state->k3);
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for y0", 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->y0);
Packit 67cb25
      free (state->k6);
Packit 67cb25
      free (state->k5);
Packit 67cb25
      free (state->k4);
Packit 67cb25
      free (state->k3);
Packit 67cb25
      free (state->k2);
Packit 67cb25
      free (state->k1);
Packit 67cb25
      free (state);
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return state;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rkck_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[], const gsl_odeiv_system * sys)
Packit 67cb25
{
Packit 67cb25
  rkck_state_t *state = (rkck_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  double *const k1 = state->k1;
Packit 67cb25
  double *const k2 = state->k2;
Packit 67cb25
  double *const k3 = state->k3;
Packit 67cb25
  double *const k4 = state->k4;
Packit 67cb25
  double *const k5 = state->k5;
Packit 67cb25
  double *const k6 = state->k6;
Packit 67cb25
  double *const ytmp = state->ytmp;
Packit 67cb25
  double *const y0 = state->y0;
Packit 67cb25
Packit 67cb25
  DBL_MEMCPY (y0, y, dim);
Packit 67cb25
Packit 67cb25
  /* k1 step */
Packit 67cb25
  if (dydt_in != NULL)
Packit 67cb25
    {
Packit 67cb25
      DBL_MEMCPY (k1, dydt_in, dim);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int s = GSL_ODEIV_FN_EVAL (sys, t, y, k1);
Packit 67cb25
      
Packit 67cb25
      if (s != GSL_SUCCESS)
Packit 67cb25
	{
Packit 67cb25
	  return s;
Packit 67cb25
	}
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    ytmp[i] = y[i] + b21 * h * k1[i];
Packit 67cb25
Packit 67cb25
  /* k2 step */
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + ah[0] * h, ytmp, k2);
Packit 67cb25
      
Packit 67cb25
    if (s != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
	return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    ytmp[i] = y[i] + h * (b3[0] * k1[i] + b3[1] * k2[i]);
Packit 67cb25
Packit 67cb25
  /* k3 step */
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + ah[1] * h, ytmp, k3);
Packit 67cb25
      
Packit 67cb25
    if (s != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
	return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    ytmp[i] = y[i] + h * (b4[0] * k1[i] + b4[1] * k2[i] + b4[2] * k3[i]);
Packit 67cb25
Packit 67cb25
  /* k4 step */
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + ah[2] * h, ytmp, k4);
Packit 67cb25
    
Packit 67cb25
    if (s != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
	return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    ytmp[i] =
Packit 67cb25
      y[i] + h * (b5[0] * k1[i] + b5[1] * k2[i] + b5[2] * k3[i] +
Packit 67cb25
                  b5[3] * k4[i]);
Packit 67cb25
Packit 67cb25
  /* k5 step */
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + ah[3] * h, ytmp, k5);
Packit 67cb25
      
Packit 67cb25
    if (s != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
	return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    ytmp[i] =
Packit 67cb25
      y[i] + h * (b6[0] * k1[i] + b6[1] * k2[i] + b6[2] * k3[i] +
Packit 67cb25
                  b6[3] * k4[i] + b6[4] * k5[i]);
Packit 67cb25
Packit 67cb25
  /* k6 step and final sum */
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + ah[4] * h, ytmp, k6);
Packit 67cb25
      
Packit 67cb25
    if (s != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
	return s;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      const double d_i = c1 * k1[i] + c3 * k3[i] + c4 * k4[i] + c6 * k6[i];
Packit 67cb25
      y[i] += h * d_i;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  /* Evaluate dydt_out[]. */
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 initial values */
Packit 67cb25
	  DBL_MEMCPY (y, y0, dim);
Packit 67cb25
	  return s;
Packit 67cb25
	}
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  /* difference between 4th and 5th order */
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      yerr[i] = h * (ec[1] * k1[i] + ec[3] * k3[i] + ec[4] * k4[i] 
Packit 67cb25
                     + ec[5] * k5[i] + ec[6] * k6[i]);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rkck_reset (void *vstate, size_t dim)
Packit 67cb25
{
Packit 67cb25
  rkck_state_t *state = (rkck_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  DBL_ZERO_MEMSET (state->k1, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k2, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k3, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k4, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k5, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->k6, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->ytmp, dim);
Packit 67cb25
  DBL_ZERO_MEMSET (state->y0, dim);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static unsigned int
Packit 67cb25
rkck_order (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rkck_state_t *state = (rkck_state_t *) vstate;
Packit 67cb25
  state = 0; /* prevent warnings about unused parameters */
Packit 67cb25
  return 5; /* FIXME: should this be 4? */
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
rkck_free (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rkck_state_t *state = (rkck_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  free (state->ytmp);
Packit 67cb25
  free (state->y0);
Packit 67cb25
  free (state->k6);
Packit 67cb25
  free (state->k5);
Packit 67cb25
  free (state->k4);
Packit 67cb25
  free (state->k3);
Packit 67cb25
  free (state->k2);
Packit 67cb25
  free (state->k1);
Packit 67cb25
  free (state);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static const gsl_odeiv_step_type rkck_type = { "rkck",  /* name */
Packit 67cb25
  1,                            /* can use dydt_in */
Packit 67cb25
  1,                            /* gives exact dydt_out */
Packit 67cb25
  &rkck_alloc,
Packit 67cb25
  &rkck_apply,
Packit 67cb25
  &rkck_reset,
Packit 67cb25
  &rkck_order,
Packit 67cb25
  &rkck_free
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
const gsl_odeiv_step_type *gsl_odeiv_step_rkck = &rkck_type;