|
Packit |
67cb25 |
/* ode-initval/rk2simp.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2004 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 |
/* Runge-Kutta 2, Gaussian implicit. Also known as implicit midpoint rule.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Non-linear equations solved by linearization, LU-decomposition
|
|
Packit |
67cb25 |
and matrix inversion. For reference, see eg.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Ascher, U.M., Petzold, L.R., Computer methods for ordinary
|
|
Packit |
67cb25 |
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_odeiv.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_linalg.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include "odeiv_util.h"
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double *Y1;
|
|
Packit |
67cb25 |
double *y0;
|
|
Packit |
67cb25 |
double *y0_orig;
|
|
Packit |
67cb25 |
double *ytmp;
|
|
Packit |
67cb25 |
double *dfdy; /* Jacobian */
|
|
Packit |
67cb25 |
double *dfdt; /* time derivatives, not used */
|
|
Packit |
67cb25 |
double *y_onestep;
|
|
Packit |
67cb25 |
gsl_permutation *p;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
rk2simp_state_t;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void *
|
|
Packit |
67cb25 |
rk2simp_alloc (size_t dim)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk2simp_state_t *state =
|
|
Packit |
67cb25 |
(rk2simp_state_t *) malloc (sizeof (rk2simp_state_t));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for rk2simp_state",
|
|
Packit |
67cb25 |
GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->Y1 = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->Y1 == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for Y1", 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->Y1);
|
|
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->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
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->ytmp = (double *) malloc (dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->ytmp == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for ytmp", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->dfdy = (double *) malloc (dim * dim * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->dfdy == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for dfdy", 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->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->dfdy);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for dfdt", 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->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
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->p = gsl_permutation_alloc (dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->p == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
free (state->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for p", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return state;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk2simp_step (double *y, rk2simp_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 2nd order semi-implicit advance with step size h.
|
|
Packit |
67cb25 |
y0 is initial values of variables y.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The linearized semi-implicit equations to calculate are:
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Y1 = y0 + h/2 * (1 - h/2 * df/dy)^(-1) * f(t + h/2, y0)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
y = y0 + h * f(t + h/2, Y1)
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const double *y0 = state->y0;
|
|
Packit |
67cb25 |
double *Y1 = state->Y1;
|
|
Packit |
67cb25 |
double *ytmp = state->ytmp;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
int s, ps;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix_view J = gsl_matrix_view_array (state->dfdy, dim, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* First solve Y1.
|
|
Packit |
67cb25 |
Calculate the inverse matrix (1 - h/2 * df/dy)^-1
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Create matrix to J */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
s = GSL_ODEIV_JA_EVAL (sys, t, y0, state->dfdy, state->dfdt);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix_scale (&J.matrix, -h / 2.0);
|
|
Packit |
67cb25 |
gsl_matrix_add_diagonal(&J.matrix, 1.0);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Invert it by LU-decomposition to invmat */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
s += gsl_linalg_LU_decomp (&J.matrix, state->p, &ps);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return GSL_EFAILED;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Evaluate f(t + h/2, y0) */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, y0, ytmp);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Calculate Y1 = y0 + h/2 * ((1-h/2 * df/dy)^-1) ytmp */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_const_view y0_view = gsl_vector_const_view_array(y0, dim);
|
|
Packit |
67cb25 |
gsl_vector_view ytmp_view = gsl_vector_view_array(ytmp, dim);
|
|
Packit |
67cb25 |
gsl_vector_view Y1_view = gsl_vector_view_array(Y1, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
s = gsl_linalg_LU_solve (&J.matrix, state->p,
|
|
Packit |
67cb25 |
&ytmp_view.vector, &Y1_view.vector);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_scale (&Y1_view.vector, 0.5 * h);
|
|
Packit |
67cb25 |
gsl_vector_add (&Y1_view.vector, &y0_view.vector);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* And finally evaluation of f(t + h/2, Y1) and calculation of y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
s = GSL_ODEIV_FN_EVAL (sys, t + 0.5 * h, Y1, ytmp);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (s != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < dim; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
y[i] = y0[i] + h * ytmp[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk2simp_apply (void *vstate, size_t dim, double t, double h,
|
|
Packit |
67cb25 |
double y[], double yerr[], const double dydt_in[],
|
|
Packit |
67cb25 |
double dydt_out[], const gsl_odeiv_system * sys)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk2simp_state_t *state = (rk2simp_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 |
/* Error estimation is done by step doubling procedure */
|
|
Packit |
67cb25 |
|
|
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 |
/* 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 = rk2simp_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 = rk2simp_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 = rk2simp_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 |
{
|
|
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]) / 3.0;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
rk2simp_reset (void *vstate, size_t dim)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk2simp_state_t *state = (rk2simp_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->Y1, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y0, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y0_orig, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->ytmp, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->dfdt, dim * dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->dfdt, dim);
|
|
Packit |
67cb25 |
DBL_ZERO_MEMSET (state->y_onestep, dim);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static unsigned int
|
|
Packit |
67cb25 |
rk2simp_order (void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk2simp_state_t *state = (rk2simp_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 |
rk2simp_free (void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
rk2simp_state_t *state = (rk2simp_state_t *) vstate;
|
|
Packit |
67cb25 |
free (state->Y1);
|
|
Packit |
67cb25 |
free (state->y0);
|
|
Packit |
67cb25 |
free (state->y0_orig);
|
|
Packit |
67cb25 |
free (state->ytmp);
|
|
Packit |
67cb25 |
free (state->dfdy);
|
|
Packit |
67cb25 |
free (state->dfdt);
|
|
Packit |
67cb25 |
free (state->y_onestep);
|
|
Packit |
67cb25 |
gsl_permutation_free (state->p);
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static const gsl_odeiv_step_type rk2simp_type = {
|
|
Packit |
67cb25 |
"rk2simp", /* name */
|
|
Packit |
67cb25 |
0, /* can use dydt_in? */
|
|
Packit |
67cb25 |
1, /* gives exact dydt_out? */
|
|
Packit |
67cb25 |
&rk2simp_alloc,
|
|
Packit |
67cb25 |
&rk2simp_apply,
|
|
Packit |
67cb25 |
&rk2simp_reset,
|
|
Packit |
67cb25 |
&rk2simp_order,
|
|
Packit |
67cb25 |
&rk2simp_free
|
|
Packit |
67cb25 |
};
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const gsl_odeiv_step_type *gsl_odeiv_step_rk2simp = &rk2simp_type;
|