Provided by: libmath-gsl-perl_0.44-1build3_amd64 

NAME
Math::GSL::ODEIV - functions for solving ordinary differential equation (ODE) initial value problems
SYNOPSIS
use Math::GSL::ODEIV qw /:all/;
DESCRIPTION
Here is a list of all the functions in this module :
• "gsl_odeiv_step_alloc($T, $dim)" - This function returns a pointer to a newly allocated instance of a
stepping function of type $T for a system of $dim dimensions.$T must be one of the step type constant
above.
• gsl_odeiv_step_reset($s) - This function resets the stepping function $s. It should be used whenever
the next use of s will not be a continuation of a previous step.
• gsl_odeiv_step_free($s) - This function frees all the memory associated with the stepping function
$s.
• gsl_odeiv_step_name($s) - This function returns a pointer to the name of the stepping function.
• gsl_odeiv_step_order($s) - This function returns the order of the stepping function on the previous
step. This order can vary if the stepping function itself is adaptive.
• "gsl_odeiv_step_apply "
• gsl_odeiv_control_alloc($T) - This function returns a pointer to a newly allocated instance of a
control function of type $T. This function is only needed for defining new types of control
functions. For most purposes the standard control functions described above should be sufficient. $T
is a gsl_odeiv_control_type.
• "gsl_odeiv_control_init($c, $eps_abs, $eps_rel, $a_y, $a_dydt) " - This function initializes the
control function c with the parameters eps_abs (absolute error), eps_rel (relative error), a_y
(scaling factor for y) and a_dydt (scaling factor for derivatives).
• "gsl_odeiv_control_free "
• "gsl_odeiv_control_hadjust "
• "gsl_odeiv_control_name "
• "gsl_odeiv_control_standard_new($eps_abs, $eps_rel, $a_y, $a_dydt)" - The standard control object is
a four parameter heuristic based on absolute and relative errors $eps_abs and $eps_rel, and scaling
factors $a_y and $a_dydt for the system state y(t) and derivatives y'(t) respectively. The step-size
adjustment procedure for this method begins by computing the desired error level D_i for each
component, D_i = eps_abs + eps_rel * (a_y |y_i| + a_dydt h |y'_i|) and comparing it with the observed
error E_i = |yerr_i|. If the observed error E exceeds the desired error level D by more than 10% for
any component then the method reduces the step-size by an appropriate factor, h_new = h_old * S *
(E/D)^(-1/q) where q is the consistency order of the method (e.g. q=4 for 4(5) embedded RK), and S is
a safety factor of 0.9. The ratio E/D is taken to be the maximum of the ratios E_i/D_i. If the
observed error E is less than 50% of the desired error level D for the maximum ratio E_i/D_i then the
algorithm takes the opportunity to increase the step-size to bring the error in line with the desired
level, h_new = h_old * S * (E/D)^(-1/(q+1)) This encompasses all the standard error scaling methods.
To avoid uncontrolled changes in the stepsize, the overall scaling factor is limited to the range 1/5
to 5.
• "gsl_odeiv_control_y_new($eps_abs, $eps_rel)" - This function creates a new control object which will
keep the local error on each step within an absolute error of $eps_abs and relative error of $eps_rel
with respect to the solution y_i(t). This is equivalent to the standard control object with a_y=1 and
a_dydt=0.
• "gsl_odeiv_control_yp_new($eps_abs, $eps_rel)" - This function creates a new control object which
will keep the local error on each step within an absolute error of $eps_abs and relative error of
$eps_rel with respect to the derivatives of the solution y'_i(t). This is equivalent to the standard
control object with a_y=0 and a_dydt=1.
• "gsl_odeiv_control_scaled_new($eps_abs, $eps_rel, $a_y, $a_dydt, $scale_abs, $dim) " - This function
creates a new control object which uses the same algorithm as gsl_odeiv_control_standard_new but with
an absolute error which is scaled for each component by the array reference $scale_abs. The formula
for D_i for this control object is, D_i = eps_abs * s_i + eps_rel * (a_y |y_i| + a_dydt h |y'_i|)
where s_i is the i-th component of the array scale_abs. The same error control heuristic is used by
the Matlab ode suite.
• gsl_odeiv_evolve_alloc($dim) - This function returns a pointer to a newly allocated instance of an
evolution function for a system of $dim dimensions.
• "gsl_odeiv_evolve_apply($e, $c, $step, $dydt, \$t, $t1, \$h, $y)" - This function advances the system
($e, $dydt) from time $t and position $y using the stepping function $step. The new time and position
are stored in $t and $y on output. The initial step-size is taken as $h, but this will be modified
using the control function $c to achieve the appropriate error bound if necessary. The routine may
make several calls to step in order to determine the optimum step-size. If the step-size has been
changed the value of $h will be modified on output. The maximum time $t1 is guaranteed not to be
exceeded by the time-step. On the final time-step the value of $t will be set to $t1 exactly.
• gsl_odeiv_evolve_reset($e) - This function resets the evolution function $e. It should be used
whenever the next use of $e will not be a continuation of a previous step.
• gsl_odeiv_evolve_free($e) - This function frees all the memory associated with the evolution function
$e.
This module also includes the following constants :
• $GSL_ODEIV_HADJ_INC
• $GSL_ODEIV_HADJ_NIL
• $GSL_ODEIV_HADJ_DEC
Step Type
• $gsl_odeiv_step_rk2 - Embedded Runge-Kutta (2, 3) method.
• $gsl_odeiv_step_rk4 - 4th order (classical) Runge-Kutta. The error estimate is obtained by halving
the step-size. For more efficient estimate of the error, use the Runge-Kutta-Fehlberg method
described below.
• $gsl_odeiv_step_rkf45 - Embedded Runge-Kutta-Fehlberg (4, 5) method. This method is a good general-
purpose integrator.
• $gsl_odeiv_step_rkck - Embedded Runge-Kutta Cash-Karp (4, 5) method.
• $gsl_odeiv_step_rk8pd - Embedded Runge-Kutta Prince-Dormand (8,9) method.
• $gsl_odeiv_step_rk2imp - Implicit 2nd order Runge-Kutta at Gaussian points.
• $gsl_odeiv_step_rk2simp
• $gsl_odeiv_step_rk4imp - Implicit 4th order Runge-Kutta at Gaussian points.
• $gsl_odeiv_step_bsimp - Implicit Bulirsch-Stoer method of Bader and Deuflhard. This algorithm
requires the Jacobian.
• $gsl_odeiv_step_gear1 - M=1 implicit Gear method.
• $gsl_odeiv_step_gear2 - M=2 implicit Gear method.
For more information on the functions, we refer you to the GSL official documentation:
<http://www.gnu.org/software/gsl/manual/html_node/>
EXAMPLE
The example is taken from <https://www.math.utah.edu/software/gsl/gsl-ref_367.html>.
use strict;
use warnings;
use Math::GSL::Errno qw($GSL_SUCCESS);
use Math::GSL::ODEIV qw/ :all /;
use Math::GSL::Matrix qw/:all/;
use Math::GSL::IEEEUtils qw/ :all /;
sub func {
my ($t, $y, $dydt, $params) = @_;
my $mu = $params->{mu};
$dydt->[0] = $y->[1];
$dydt->[1] = -$y->[0] - $mu*$y->[1]*(($y->[0])**2 - 1);
return $GSL_SUCCESS;
}
sub jac {
my ($t, $y, $dfdy, $dfdt, $params) = @_;
my $mu = $params->{mu};
my $m = gsl_matrix_view_array($dfdy, 2, 2);
gsl_matrix_set( $m, 0, 0, 0.0 );
gsl_matrix_set( $m, 0, 1, 1.0 );
gsl_matrix_set( $m, 1, 0, (-2.0 * $mu * $y->[0] * $y->[1]) - 1.0 );
gsl_matrix_set( $m, 1, 1, -$mu * (($y->[0])**2 - 1.0) );
$dfdt->[0] = 0.0;
$dfdt->[1] = 0.0;
return $GSL_SUCCESS;
}
my $T = $gsl_odeiv_step_rk8pd;
my $s = gsl_odeiv_step_alloc($T, 2);
my $c = gsl_odeiv_control_y_new(1e-6, 0.0);
my $e = gsl_odeiv_evolve_alloc(2);
my $params = { mu => 10 };
my $sys = Math::GSL::ODEIV::gsl_odeiv_system->new(\&func, \&jac, 2, $params );
my $t = 0.0;
my $t1 = 100.0;
my $h = 1e-6;
my $y = [ 1.0, 0.0 ];
gsl_ieee_env_setup;
while ($t < $t1) {
my $status = gsl_odeiv_evolve_apply ($e, $c, $s, $sys, \$t, $t1, \$h, $y);
last if $status != $GSL_SUCCESS;
printf "%.5e %.5e %.5e\n", $t, $y->[0], $y->[1];
}
gsl_odeiv_evolve_free($e);
gsl_odeiv_control_free($c);
gsl_odeiv_step_free($s);
AUTHORS
Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
COPYRIGHT AND LICENSE
Copyright (C) 2008-2023 Jonathan "Duke" Leto and Thierry Moisan
This program is free software; you can redistribute it and/or modify it under the same terms as Perl
itself.
perl v5.38.2 2024-03-31 Math::GSL::ODEIV(3pm)