Blame ode-initval2/rk2.c

Packit 67cb25
/* ode-initval/rk2.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 2(3), Euler-Cauchy */
Packit 67cb25
Packit 67cb25
/* Author:  G. Jungman
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
/* Reference: Abramowitz & Stegun, section 25.5. Runge-Kutta 2nd (25.5.7)
Packit 67cb25
   and 3rd (25.5.8) order methods */
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_odeiv2.h>
Packit 67cb25
Packit 67cb25
#include "odeiv_util.h"
Packit 67cb25
#include "step_utils.c"
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  double *k1;
Packit 67cb25
  double *k2;
Packit 67cb25
  double *k3;
Packit 67cb25
  double *ytmp;
Packit 67cb25
}
Packit 67cb25
rk2_state_t;
Packit 67cb25
Packit 67cb25
static void *
Packit 67cb25
rk2_alloc (size_t dim)
Packit 67cb25
{
Packit 67cb25
  rk2_state_t *state = (rk2_state_t *) malloc (sizeof (rk2_state_t));
Packit 67cb25
Packit 67cb25
  if (state == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for rk2_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->ytmp = (double *) malloc (dim * sizeof (double));
Packit 67cb25
Packit 67cb25
  if (state->ytmp == 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 ytmp", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return state;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rk2_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_odeiv2_system * sys)
Packit 67cb25
{
Packit 67cb25
  rk2_state_t *state = (rk2_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 ytmp = state->ytmp;
Packit 67cb25
Packit 67cb25
  /* k1 step */
Packit 67cb25
  /* k1 = f(t,y) */
Packit 67cb25
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
  /* k2 step */
Packit 67cb25
  /* k2 = f(t + 0.5*h, y + 0.5*k1) */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      ytmp[i] = y[i] + 0.5 * h * k1[i];
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * 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
  /* k3 step */
Packit 67cb25
  /* for 3rd order estimates, is used for error estimation
Packit 67cb25
     k3 = f(t + h, y - k1 + 2*k2) */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      ytmp[i] = y[i] + h * (-k1[i] + 2.0 * k2[i]);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int s = GSL_ODEIV_FN_EVAL (sys, t + 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
  /* final sum */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < dim; i++)
Packit 67cb25
    {
Packit 67cb25
      /* Save original values if derivative evaluation below fails */
Packit 67cb25
      ytmp[i] = y[i];
Packit 67cb25
Packit 67cb25
      {
Packit 67cb25
        const double ksum3 = (k1[i] + 4.0 * k2[i] + k3[i]) / 6.0;
Packit 67cb25
        y[i] += h * ksum3;
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, ytmp, dim);
Packit 67cb25
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
      const double ksum3 = (k1[i] + 4.0 * k2[i] + k3[i]) / 6.0;
Packit 67cb25
      yerr[i] = h * (k2[i] - ksum3);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
rk2_reset (void *vstate, size_t dim)
Packit 67cb25
{
Packit 67cb25
  rk2_state_t *state = (rk2_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->ytmp, dim);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static unsigned int
Packit 67cb25
rk2_order (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rk2_state_t *state = (rk2_state_t *) vstate;
Packit 67cb25
  state = 0;                    /* prevent warnings about unused parameters */
Packit 67cb25
  return 2;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
rk2_free (void *vstate)
Packit 67cb25
{
Packit 67cb25
  rk2_state_t *state = (rk2_state_t *) vstate;
Packit 67cb25
  free (state->k1);
Packit 67cb25
  free (state->k2);
Packit 67cb25
  free (state->k3);
Packit 67cb25
  free (state->ytmp);
Packit 67cb25
  free (state);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static const gsl_odeiv2_step_type rk2_type = { "rk2",   /* name */
Packit 67cb25
  1,                            /* can use dydt_in */
Packit 67cb25
  1,                            /* gives exact dydt_out */
Packit 67cb25
  &rk2_alloc,
Packit 67cb25
  &rk2_apply,
Packit 67cb25
  &stepper_set_driver_null,
Packit 67cb25
  &rk2_reset,
Packit 67cb25
  &rk2_order,
Packit 67cb25
  &rk2_free
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
const gsl_odeiv2_step_type *gsl_odeiv2_step_rk2 = &rk2_type;