view scripts/ode/private/runge_kutta_23.m @ 22626:869c02fde46c stable

Further clean-up of ode functions. * scripts/ode/private/AbsRel_Norm.m: Renamed to AbsRel_norm.m * scripts/ode/module.mk: Add AbsRel_norm.m to build system. * ode23.m, ode45.m: Remove extra comma from Copyright statement. Add ode45 to @seealso links in docstring. Add FIXME notes for questionable code. Use numel instead of length for clarity. Use space after function name and opening parenthesis. Wrap long lines to less than 80 characters. Change odemergeopts function call to match new prototype. Use single quotes to simplify strings that contain double quotes. Use 2-space indent in %!function blocks. * odeget.m: Remove extra comma from Copyright statement. Use double quote in preference to single quote. Delete whitespace at end of lines. * odeplot.m: Re-write docstring. Include @seealso links in docstring. Declare all persistent variables in a single declaration. Use in-place += operator for efficiency. Add FIXME notes for questionable code. * odeset.m: Remove extra comma from Copyright statement. Add additional calling form with 1 output and 0 inputs to docstring. Add FIXME notes for questionable code. Delete whitespace at end of lines. * odeset.m(print_options): Use single quotes to simplify strings with double quotes. Put default value of option first in list. * integrate_adaptive.m: Wrap long lines < 80 characters. Delete whitespace at end of lines. Correct indentation of declared values after '='. * kahan.m: Reise docstring. * ode_event_handler.m: Use retval in docstring to match functin prototype. Delete unnecessary comments. * odedefaults.m: Remove extra comma from Copyright statement. Match variable names in docstring to function prototype. Use space after function name and before opening parenthesis. Delete whitespace at end of lines. * odemergeopts.m: Remove extra comma from Copyright statement. Delete extra space in function prototype. Change function prototype to have "caller" as first argument to match rest of Octave. * runge_kutta_23.m: Clean up declaration of persistent variables. Correct misspellings in comments. * runge_kutta_45_dorpri.m: Put description of input arguments before output arguments in docstring. Clean up declaration of persistent variables. * runge_kutta_interpolate.m: Use double quotes in preference to single quotes. Eliminate line continuations for code that would fit on a single line. Remove obsolete code that calls non-existent functions. Capitalize Hermite in comments. Cleanup declaration of persistent variables. * starting_stepsize.m: Replace calls to AbsRel_Norm with AbsRel_norm.
author Rik <rik@octave.org>
date Sat, 15 Oct 2016 22:39:29 -0700
parents bac0d6f07a3e
children 3a2b891d0b33 e9a0469dedd9
line wrap: on
line source

## Copyright (C) 2013-2016 Jacopo Corno <jacopo.corno@gmail.com>
## Copyright (C) 2013-2016 Roberto Porcu' <roberto.porcu@polimi.it>
##
## 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  {} {[@var{t_next}, @var{x_next}] =} runge_kutta_23 (@var{fun}, @var{t}, @var{x}, @var{dt})
## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_23 (@var{fun}, @var{t}, @var{x}, @var{dt}, @var{options})
## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_23 (@var{fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals_in})
## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_23 (@var{fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals_in}, @var{t_next})
## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}] =} runge_kutta_23 (@dots{})
## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}, @var{k_vals_out}] =} runge_kutta_23 (@dots{})
##
## 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 Bogacki-Shampine
## method of third order.  For the definition of this method see
## @url{http://en.wikipedia.org/wiki/List_of_Runge%E2%80%93Kutta_methods}.
##
## @var{fun} is a function handle that defines the ODE: @code{y' = f(tau,y)}.
## The function must accept two inputs where the first is time @var{tau} and
## the second is a column vector of unknowns @var{y}.
##
## @var{t} is the first extreme of integration interval.
##
## @var{x} is the initial condition of the system..
##
## @var{dt} is the timestep, that is the length of the integration interval.
##
## The optional fourth argument @var{options} specifies options for the ODE
## solver.  It is a structure generated by @code{odeset}.  In particular it
## contains the field @var{funarguments} with the optional arguments to be used
## in the evaluation of @var{fun}.
##
## The optional fifth argument @var{k_vals_in} contains the Runge-Kutta
## evaluations of the previous step to use in a FSAL scheme.
##
## The optional sixth argument @var{t_next} (@code{t_next = t + dt}) specifies
## the end of the integration interval.  The output @var{x_next} s the higher
## order computed solution at time @var{t_next} (local extrapolation is
## performed).
##
## Optionally the functions can also return @var{x_est}, a lower order solution
## for the estimation of the error, and @var{k_vals_out}, a matrix containing
## the Runge-Kutta evaluations to use in a FSAL scheme or for dense output.
##
## @seealso{runge_kutta_45_dorpri}
## @end deftypefn

function [t_next, x_next, x_est, k] = runge_kutta_23 (f, t, x, dt,
                                                      options = [],
                                                      k_vals = [],
                                                      t_next = t + dt)

  persistent a = [0   0   0;
                  1/2 0   0;
                  0   3/4 0];
  persistent b = [0, 1/2, 3/4, 1];
  persistent c = [2/9, 1/3, 4/9];
  persistent c_prime = [7/24, 1/4, 1/3, 1/8];

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

  if (! isempty (options))  # extra arguments for function evaluator
    args = options.funarguments;
  else
    args = {};
  endif

  if (! isempty (k_vals))    # k values from previous step are passed
    k(:,1) = k_vals(:,end);  # FSAL property
  else
    k(:,1) = feval (f, t, x, args{:});
  endif

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

  ## compute new time and new values for the unknowns
  ## t_next = t + dt;
  x_next = x + k(:,1:3) * cc(:);  # 3rd order approximation

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

endfunction