Mercurial > octave
diff scripts/ode/ode15i.m @ 25026:f886561f9696 stable
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
* NEWS: Announce new functions decic, ode15i, ode15s.
* diffeq.txi: Add decic, ode15i, ode15s to manual. Rewrite table of solvers.
* decic.m, ode15i.m, ode15s.m, ode23.m, ode45.m: Rewrite documentation.
author | Colin Macdonald <cbm@m.fsf.org> |
---|---|
date | Tue, 27 Mar 2018 00:52:29 -0700 |
parents | 194eb4bd202b |
children | 6652d3823428 |
line wrap: on
line diff
--- a/scripts/ode/ode15i.m Tue Mar 27 14:03:52 2018 -0700 +++ b/scripts/ode/ode15i.m Tue Mar 27 00:52:29 2018 -0700 @@ -23,16 +23,16 @@ ## @deftypefnx {} {@var{solution} =} ode15i (@dots{}) ## @deftypefnx {} {} ode15i (@dots{}) ## -## Solve a set of full-implicit Ordinary Differential Equations and -## Differential Algebraic Equations (DAEs) of index 1, with the variable-step, -## variable order BDF (Backward Differentiation Formula) method, which -## ranges from order 1 to 5. +## Solve a set of fully-implicit Ordinary Differential Equations (ODEs) or +## Differential Algebraic Equations (DAEs) of index 1, with a variable step, +## variable order BDF (Backward Differentiation Formula) method that ranges +## from order 1 to 5. ## ## @var{fun} is a function handle, inline function, or string containing the -## name of the function that defines the ODE: @code{y' = f(t,y,yp)}. The +## name of the function that defines the ODE: @code{0 = f(t,y,yp)}. The ## function must accept three inputs where the first is time @var{t}, the -## second is a column vector of unknowns @var{y}, and the third is a column -## vector of unknowns @var{yp}. +## second is the function value @var{y} (a column vector), and the third +## is the derivative value @var{yp} (a column vector). ## ## @var{trange} specifies the time interval over which the ODE will be ## evaluated. Typically, it is a two-element vector specifying the initial and @@ -46,8 +46,8 @@ ## value in @var{y0} and @var{yp0}. ## ## @var{y0} and @var{yp0} must be consistent initial conditions, meaning that -## @code{f(t,y0,yp0) = 0} is satisfied. You can use function @code{decic} to -## compute consistent initial conditions, given initial guesses. +## @code{f(t,y0,yp0) = 0} is satisfied. The function @code{decic} may be used +## to compute consistent initial conditions given initial guesses. ## ## The optional fifth argument @var{ode_opt} specifies non-default options to ## the ODE solver. It is a structure generated by @code{odeset}. @@ -57,46 +57,46 @@ ## output @var{y} is a matrix in which each column refers to a different ## unknown of the problem and each row corresponds to a time in @var{t}. ## -## The output can also be returned as a structure @var{solution} which -## has a field @var{x} containing a row vector of times where the solution -## was evaluated and a field @var{y} containing the solution matrix such -## that each column corresponds to a time in @var{x}. -## Use @code{fieldnames (@var{solution})} to see the other fields and +## The output can also be returned as a structure @var{solution} which has a +## field @var{x} containing a row vector of times where the solution was +## evaluated and a field @var{y} containing the solution matrix such that each +## column corresponds to a time in @var{x}. Use +## @w{@code{fieldnames (@var{solution})}} to see the other fields and ## additional information returned. ## -## If no output arguments are requested, and no @code{OutputFcn} is -## specified in @var{ode_opt}, then the @code{OutputFcn} is set to -## @code{odeplot} and the results of the solver are plotted immediately. +## If no output arguments are requested, and no @code{OutputFcn} is specified +## in @var{ode_opt}, then the @code{OutputFcn} is set to @code{odeplot} and the +## results of the solver are plotted immediately. ## -## If using the @qcode{"Events"} option then three additional outputs may -## be returned. @var{te} holds the time when an Event function returned a -## zero. @var{ye} holds the value of the solution at time @var{te}. @var{ie} +## If using the @qcode{"Events"} option then three additional outputs may be +## returned. @var{te} holds the time when an Event function returned a zero. +## @var{ye} holds the value of the solution at time @var{te}. @var{ie} ## contains an index indicating which Event function was triggered in the case ## of multiple Event functions. ## -## Example: Solve the @nospell{Robetson's} equations: +## Example: Solve @nospell{Robertson's} equations: ## -## @example +## @smallexample ## @group -## function r = robertsidae (@var{t}, @var{y}, @var{yp}) -## r = [-(@var{yp}(1) + 0.04*@var{y}(1) - 1e4*@var{y}(2)*@var{y}(3)); -## -(@var{yp}(2) - 0.04*@var{y}(1) + 1e4*@var{y}(2)*@var{y}(3) + 3e7*@var{y}(2)^2); -## @var{y}(1) + @var{y}(2) + @var{y}(3) - 1]; +## function r = robertson_dae (@var{t}, @var{y}, @var{yp}) +## r = [ -(@var{yp}(1) + 0.04*@var{y}(1) - 1e4*@var{y}(2)*@var{y}(3)) +## -(@var{yp}(2) - 0.04*@var{y}(1) + 1e4*@var{y}(2)*@var{y}(3) + 3e7*@var{y}(2)^2) +## @var{y}(1) + @var{y}(2) + @var{y}(3) - 1 ]; ## endfunction -## [@var{t},@var{y}] = ode15i (@@robertsidae, [0, 1e3], [1; 0; 0], [-1e-4; 1e-4; 0]); +## [@var{t},@var{y}] = ode15i (@@robertson_dae, [0, 1e3], [1; 0; 0], [-1e-4; 1e-4; 0]); ## @end group -## @end example +## @end smallexample ## @seealso{decic, odeset, odeget} ## @end deftypefn function varargout = ode15i (fun, trange, y0, yp0, varargin) - solver = "ode15i"; - if (nargin < 4) print_usage (); endif + solver = "ode15i"; + n = numel (y0); if (nargin > 4)