Math::GSL::ODEIV(3pm) User Contributed Perl Documentation Math::GSL::ODEIV(3pm)

Math::GSL::ODEIV - functions for solving ordinary differential equation (ODE) initial value problems

 use Math::GSL::ODEIV qw /:all/;

Here is a list of all the functions in this module :

This module also includes the following constants :

  • $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/>

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);

Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>

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.

2024-03-31 perl v5.38.2