view scripts/ode/private/runge_kutta_45_dorpri.m @ 20584:eb9e2d187ed2

maint: Use Octave coding conventions in scripts/ode/private dir. * AbsRel_Norm.m, fuzzy_compare.m, hermite_quartic_interpolation.m, integrate_adaptive.m, integrate_const.m, integrate_n_steps.m, kahan.m, ode_struct_value_check.m, odepkg_event_handle.m, odepkg_structure_check.m, runge_kutta_45_dorpri.m, starting_stepsize.m: Wrap long lines to < 80 chars. Use double quotes rather than single quotes where possible. Use ';' at end of keywords "return;" and "break;" Use '##" for stand-alone comments and '#' for end-of-line comments. Use two spaces after period before starting new sentence. Use '!' instead of '~' for logical negation. Use specific form of end (endif, endfor, etc.). Don't use line continuation marker '...' unless necessary.
author Rik <rik@octave.org>
date Sun, 04 Oct 2015 22:18:54 -0700
parents 25623ef2ff4f
children b7ac1e94266e
line wrap: on
line source

## Copyright (C) 2013, Roberto Porcu' <roberto.porcu@polimi.it>
## Copyright (C) 2015, Carlo de Falco
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation; either version 3 of the License, or (at
## your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
## General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals_in})
## @deftypefnx {Function File} {[@var{t_next}, @var{x_next}, @var{x_est}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals_in})
## @deftypefnx {Function File} {[@var{t_next}, @var{x_next}, @var{x_est}, @var{k_vals_out}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals_in})
##
## This function can be used to integrate a system of ODEs with a given initial
## condition @var{x} from @var{t} to @var{t+dt}, with the Dormand-Prince method.
## For the definition of this method see
## @url{http://en.wikipedia.org/wiki/Dormand%E2%80%93Prince_method}.
##
## First output argument is the final integration time value.
##
## Second output parameter is the higher order computed solution at time
## @var{t_next} (local extrapolation).
##
## Third output parameter is a lower order solution for the estimation of the
## error.
##
## Fourth output parameter is matrix containing the Runge-Kutta evaluations
## to use in a FSAL scheme or for dense output.
##
## First input argument is the function describing the system of ODEs to be
## integrated.
##
## Second input parameter is the first extreme of integration interval.
##
## Third input argument is the initial condition of the system.
##
## Fourth input argument is the timestep, that is the length of the
## integration interval.
##
## Fifth input parameter is optional and describes a set of options useful to
## adapt the computation to what is needed.
##
## Sixth input parameter is optional and describes the Runge-Kutta evaluations
## of the previous step to use in a FSAL scheme.
## @end deftypefn
##
## @seealso{odepkg}

function [t_out, x_out, x_est, k] = ...
         runge_kutta_45_dorpri (f, t, x, dt, varargin)

  persistent a = [0           0          0           0        0          0;
                  1/5         0          0           0        0          0;
                  3/40        9/40       0           0        0          0;
                  44/45      -56/15      32/9        0        0          0;
                  19372/6561 -25360/2187 64448/6561 -212/729  0          0;
                  9017/3168  -355/33     46732/5247  49/176  -5103/18656 0;
                  35/384      0          500/1113    125/192 -2187/6784  11/84];
  persistent b = [0 1/5 3/10 4/5 8/9 1 1];
  persistent c = [(35/384) 0 (500/1113) (125/192) (-2187/6784) (11/84)];
  ## x_est according to Shampine 1986:
  ## persistent c_prime = [(1951/21600) 0 (22642/50085) (451/720), ...
  ##                       (-12231/42400) (649/6300) (1/60)];
  persistent c_prime = [(5179/57600) 0 (7571/16695) (393/640), ...
                        (-92097/339200) (187/2100)  (1/40)];

  s = t + dt * b;
  cc = dt * c;
  aa = dt * a;
  k = zeros (rows (x), 7);

  if (nargin >= 5) # options are passed
    args = varargin{1}.vfunarguments;
    if (nargin >= 6) # both the options and the k values are passed
      k(:,1) = varargin{2}(:,end); # FSAL property
    else      
      k(:,1) = feval (f, t, x, args{:});
    endif
  else
    args = {};
  endif

  k(:,2) = feval (f, s(2), x + k(:,1) * aa(2, 1).', args{:});
  k(:,3) = feval (f, s(3), x + k(:,1:2) * aa(3, 1:2).', args{:});
  k(:,4) = feval (f, s(4), x + k(:,1:3) * aa(4, 1:3).', args{:});
  k(:,5) = feval (f, s(5), x + k(:,1:4) * aa(5, 1:4).', args{:});
  k(:,6) = feval (f, s(6), x + k(:,1:5) * aa(6, 1:5).', args{:});

  ## compute new time and new values for the unknowns
  t_out = t + dt;
  x_out = x + k(:,1:6) * cc(:);  # 5th order approximation

  ## if the estimation of the error is required
  if (nargout >= 3)
    ## new solution to be compared with the previous one
    k(:,7) = feval (f, t + dt, x_out, args{:});
    cc_prime = dt * c_prime;
    x_est = x + k * cc_prime(:);
  endif
  
endfunction