annotate scripts/ode/ode15s.m @ 28041:5e44268dca6f

Replace input validation relying on str2func with alternatives (bug #57351). * ode15i.m, ode15s.m, ode23.m, ode45.m, check_default_input.m: Replace try/catch blocks around str2func with a call to exist() to check whether function exists before calling str2func.
author Rik <rik@octave.org>
date Mon, 03 Feb 2020 20:22:19 -0800
parents bd51beb6205e
children ace8dc642e4d
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
1 ########################################################################
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
2 ##
27919
1891570abac8 update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 27918
diff changeset
3 ## Copyright (C) 2016-2020 The Octave Project Developers
27918
b442ec6dda5c use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
4 ##
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
5 ## See the file COPYRIGHT.md in the top-level directory of this
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
6 ## distribution or <https://octave.org/copyright/>.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
7 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
8 ## This file is part of Octave.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
9 ##
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24012
diff changeset
10 ## Octave is free software: you can redistribute it and/or modify it
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
11 ## under the terms of the GNU General Public License as published by
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24012
diff changeset
12 ## the Free Software Foundation, either version 3 of the License, or
24007
e8a74d95b4f3 maint: Use same format for Copyright statement throught code base.
Rik <rik@octave.org>
parents: 23023
diff changeset
13 ## (at your option) any later version.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
14 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
15 ## Octave is distributed in the hope that it will be useful, but
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of
24007
e8a74d95b4f3 maint: Use same format for Copyright statement throught code base.
Rik <rik@octave.org>
parents: 23023
diff changeset
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
e8a74d95b4f3 maint: Use same format for Copyright statement throught code base.
Rik <rik@octave.org>
parents: 23023
diff changeset
18 ## GNU General Public License for more details.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
19 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
20 ## You should have received a copy of the GNU General Public License
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
21 ## along with Octave; see the file COPYING. If not, see
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24012
diff changeset
22 ## <https://www.gnu.org/licenses/>.
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
23 ##
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
24 ########################################################################
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
25
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
26 ## -*- texinfo -*-
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
27 ## @deftypefn {} {[@var{t}, @var{y}] =} ode15s (@var{fun}, @var{trange}, @var{y0})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
28 ## @deftypefnx {} {[@var{t}, @var{y}] =} ode15s (@var{fun}, @var{trange}, @var{y0}, @var{ode_opt})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
29 ## @deftypefnx {} {[@var{t}, @var{y}, @var{te}, @var{ye}, @var{ie}] =} ode15s (@dots{})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
30 ## @deftypefnx {} {@var{solution} =} ode15s (@dots{})
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
31 ## @deftypefnx {} {} ode15s (@dots{})
25113
476fc012d09a doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents: 25054
diff changeset
32 ## Solve a set of stiff Ordinary Differential Equations (ODEs) or stiff
476fc012d09a doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents: 25054
diff changeset
33 ## semi-explicit index 1 Differential Algebraic Equations (DAEs).
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
34 ##
25113
476fc012d09a doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents: 25054
diff changeset
35 ## @code{ode15s} uses a variable step, variable order BDF (Backward
476fc012d09a doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents: 25054
diff changeset
36 ## Differentiation Formula) method that ranges from order 1 to 5.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
37 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
38 ## @var{fun} is a function handle, inline function, or string containing the
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
39 ## name of the function that defines the ODE: @code{y' = f(t,y)}. The function
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
40 ## must accept two inputs where the first is time @var{t} and the second is a
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
41 ## column vector of unknowns @var{y}.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
42 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
43 ## @var{trange} specifies the time interval over which the ODE will be
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
44 ## evaluated. Typically, it is a two-element vector specifying the initial and
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
45 ## final times (@code{[tinit, tfinal]}). If there are more than two elements
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
46 ## then the solution will also be evaluated at these intermediate time
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
47 ## instances.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
48 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
49 ## @var{init} contains the initial value for the unknowns. If it is a row
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
50 ## vector then the solution @var{y} will be a matrix in which each column is
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
51 ## the solution for the corresponding initial value in @var{init}.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
52 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
53 ## The optional fourth argument @var{ode_opt} specifies non-default options to
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
54 ## the ODE solver. It is a structure generated by @code{odeset}.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
55 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
56 ## The function typically returns two outputs. Variable @var{t} is a
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
57 ## column vector and contains the times where the solution was found. The
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
58 ## output @var{y} is a matrix in which each column refers to a different
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
59 ## unknown of the problem and each row corresponds to a time in @var{t}.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
60 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
61 ## The output can also be returned as a structure @var{solution} which has a
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
62 ## field @var{x} containing a row vector of times where the solution was
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
63 ## evaluated and a field @var{y} containing the solution matrix such that each
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
64 ## column corresponds to a time in @var{x}. Use
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
65 ## @w{@code{fieldnames (@var{solution})}} to see the other fields and
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
66 ## additional information returned.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
67 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
68 ## If no output arguments are requested, and no @code{OutputFcn} is specified
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
69 ## in @var{ode_opt}, then the @code{OutputFcn} is set to @code{odeplot} and the
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
70 ## results of the solver are plotted immediately.
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
71 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
72 ## If using the @qcode{"Events"} option then three additional outputs may be
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
73 ## returned. @var{te} holds the time when an Event function returned a zero.
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
74 ## @var{ye} holds the value of the solution at time @var{te}. @var{ie}
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
75 ## contains an index indicating which Event function was triggered in the case
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
76 ## of multiple Event functions.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
77 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
78 ## Example: Solve @nospell{Robertson's} equations:
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
79 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
80 ## @smallexample
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
81 ## @group
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
82 ## function r = robertson_dae (@var{t}, @var{y})
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
83 ## r = [ -0.04*@var{y}(1) + 1e4*@var{y}(2)*@var{y}(3)
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
84 ## 0.04*@var{y}(1) - 1e4*@var{y}(2)*@var{y}(3) - 3e7*@var{y}(2)^2
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
85 ## @var{y}(1) + @var{y}(2) + @var{y}(3) - 1 ];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
86 ## endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
87 ## opt = odeset ("Mass", [1 0 0; 0 1 0; 0 0 0], "MStateDependence", "none");
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
88 ## [@var{t},@var{y}] = ode15s (@@robertson_dae, [0, 1e3], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
89 ## @end group
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
90 ## @end smallexample
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
91 ## @seealso{decic, odeset, odeget, ode23, ode45}
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
92 ## @end deftypefn
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
93
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
94 function varargout = ode15s (fun, trange, y0, varargin)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
95
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
96 if (nargin < 3)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
97 print_usage ();
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
98 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
99
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
100 solver = "ode15s";
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
101 ## Check fun, trange, y0, yp0
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
102 fun = check_default_input (fun, trange, solver, y0);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
103
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
104 n = numel (y0);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
105
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
106 if (nargin > 3)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
107 options = varargin{1};
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
108 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
109 options = odeset ();
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
110 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
111
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
112 if (! isempty (options.Mass))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
113 if (ischar (options.Mass))
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
114 if (! exist (options.Mass))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
115 error ("Octave:invalid-input-arg",
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
116 [solver ": function '" options.Mass "' not found"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
117 endif
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
118 options.Mass = str2func (options.Mass);
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
119 endif
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
120 if (! is_function_handle (options.Mass))
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
121 error ("Octave:invalid-input-arg",
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
122 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
123 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
124 endif
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
125
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
126 if (! isempty (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
127 if (ischar (options.Jacobian))
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
128 if (! exist (options.Jacobian))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
129 error ("Octave:invalid-input-arg",
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
130 [solver ": function '" options.Jacobian "' not found"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
131 endif
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
132 options.Jacobian = str2func (options.Jacobian);
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
133 endif
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
134 if (! is_function_handle (options.Jacobian))
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
135 error ("Octave:invalid-input-arg",
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
136 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
137 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
138 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
139
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
140 if (! isempty (options.OutputFcn))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
141 if (ischar (options.OutputFcn))
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
142 if (! exist (options.OutputFcn))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
143 error ("Octave:invalid-input-arg",
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
144 [solver ": function '" options.OutputFcn "' not found"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
145 endif
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
146 options.OutputFcn = str2func (options.OutputFcn);
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
147 endif
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
148 if (! is_function_handle (options.OutputFcn))
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
149 error ("Octave:invalid-input-arg",
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
150 [solver ": invalid value assigned to field 'OutputFcn'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
151 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
152 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
153
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
154 if (! isempty (options.Events))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
155 if (ischar (options.Events))
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
156 if (! exist (options.Events))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
157 error ("Octave:invalid-input-arg",
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
158 [solver ": function '" options.Events "' not found"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
159 endif
28041
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
160 options.Events = str2func (options.Events);
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
161 endif
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
162 if (! is_function_handle (options.Events))
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
163 error ("Octave:invalid-input-arg",
5e44268dca6f Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents: 27923
diff changeset
164 [solver ": invalid value assigned to field 'Events'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
165 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
166 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
167
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
168 [defaults, classes, attributes] = odedefaults (n, trange(1), trange(end));
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
169
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
170 classes = odeset (classes, "Vectorized", {});
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
171 attributes = odeset (attributes, "Jacobian", {}, "Vectorized", {});
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
172
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
173 options = odemergeopts ("ode15s", options, defaults,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
174 classes, attributes, solver);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
175
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
176 ## Mass
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
177 options.havemassfun = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
178 options.havestatedep = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
179 options.havetimedep = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
180 options.havemasssparse = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
181
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
182 if (! isempty (options.Mass))
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25113
diff changeset
183 if (is_function_handle (options.Mass))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
184 options.havemassfun = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
185 if (nargin (options.Mass) == 2)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
186 options.havestatedep = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
187 M = options.Mass (trange(1), y0);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
188 options.havemasssparse = issparse (M);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
189 if (any (size (M) != [n n]) || ! isnumeric (M) || ! isreal (M))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
190 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
191 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
192 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
193 elseif (nargin (options.Mass) == 1)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
194 options.havetimedep = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
195 M = options.Mass (trange(1));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
196 options.havemasssparse = issparse (M);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
197 if (any (size (M) != [n n]) || ! isnumeric (M) || ! isreal (M))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
198 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
199 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
200 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
201 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
202 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
203 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
204 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
205 elseif (ismatrix (options.Mass))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
206 options.havemasssparse = issparse (options.Mass);
25949
c2bf210ac94f style fixes
John W. Eaton <jwe@octave.org>
parents: 25113
diff changeset
207 if (any (size (options.Mass) != [n n])
c2bf210ac94f style fixes
John W. Eaton <jwe@octave.org>
parents: 25113
diff changeset
208 || ! isnumeric (options.Mass) || ! isreal (options.Mass))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
209 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
210 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
211 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
212 else
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
213 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
214 [solver ": invalid value assigned to field 'Mass'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
215 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
216 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
217
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
218 ## Jacobian
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
219 options.havejac = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
220 options.havejacsparse = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
221 options.havejacfun = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
222
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
223 if (! isempty (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
224 options.havejac = true;
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25113
diff changeset
225 if (is_function_handle (options.Jacobian))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
226 options.havejacfun = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
227 if (nargin (options.Jacobian) == 2)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
228 [A] = options.Jacobian (trange(1), y0);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
229 if (issparse (A))
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
230 options.havejacsparse = true; # Jac is sparse fun
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
231 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
232 if (any (size (A) != [n n]) || ! isnumeric (A) || ! isreal (A))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
233 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
234 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
235 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
236 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
237 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
238 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
239 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
240 elseif (ismatrix (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
241 if (issparse (options.Jacobian))
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
242 options.havejacsparse = true; # Jac is sparse matrix
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
243 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
244 if (! issquare (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
245 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
246 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
247 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
248 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
249 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
250 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
251 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
252 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
253
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
254 ## Derivative of M(t,y) for implicit problem not implemented yet
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
255 if (! isempty (options.Mass) && ! isempty (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
256 if (options.MStateDependence != "none" || options.havestatedep == true)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
257 options.havejac = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
258 options.Jacobian = [];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
259 warning ("ode15s:mass_state_dependent_provided",
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
260 ["with MStateDependence != 'none' an internal", ...
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
261 " approximation of Jacobian Matrix will be used.", ...
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
262 " Set MStateDependence equal to 'none' if you want ", ...
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
263 " to provide a constant or time-dependent Jacobian"]);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
264 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
265 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
266
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
267 ## Use sparse methods only if all matrices are sparse
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
268 if (! options.havemasssparse)
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
269 options.havejacsparse = false;
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
270 endif
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
271
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
272 ## If Mass or Jacobian is fun, then new Jacobian is fun
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
273 if (options.havejac)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
274 if (options.havejacfun || options.havetimedep)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
275 options.Jacobian = @ (t, y, yp) wrapjacfun (t, y, yp,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
276 options.Jacobian,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
277 options.Mass,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
278 options.havetimedep,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
279 options.havejacfun);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
280 options.havejacfun = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
281 else ## All matrices are constant
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
282 options.Jacobian = {[- options.Jacobian], [options.Mass]};
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
283
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
284 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
285 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
286
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
287 ## Abstol and Reltol
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
288 options.haveabstolvec = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
289
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
290 if (numel (options.AbsTol) != 1 && numel (options.AbsTol) != n)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
291 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
292 [solver ": invalid value assigned to field 'AbsTol'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
293 elseif (numel (options.AbsTol) == n)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
294 options.haveabstolvec = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
295 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
296
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
297 ## Stats
22935
c9344df03da5 Allow case-insensitive option argument 'on' to ode solvers (bug #49918).
Rik <rik@octave.org>
parents: 22933
diff changeset
298 options.havestats = strcmpi (options.Stats, "on");
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
299
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
300 ## Don't use Refine when the output is a structure
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
301 if (nargout == 1)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
302 options.Refine = 1;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
303 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
304
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
305 ## OutputFcn and OutputSel
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
306 if (isempty (options.OutputFcn) && nargout == 0)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
307 options.OutputFcn = @odeplot;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
308 options.haveoutputfunction = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
309 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
310 options.haveoutputfunction = ! isempty (options.OutputFcn);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
311 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
312
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
313 options.haveoutputselection = ! isempty (options.OutputSel);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
314 if (options.haveoutputselection)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
315 options.OutputSel = options.OutputSel - 1;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
316 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
317
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
318 ## Events
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
319 options.haveeventfunction = ! isempty (options.Events);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
320
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
321 yp0 = options.InitialSlope;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
322
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
323 [t, y, te, ye, ie] = __ode15__ (@ (t, y, yp) wrap (t, y, yp, options.Mass,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
324 options.havetimedep,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
325 options.havestatedep,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
326 fun),
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
327 trange, y0, yp0, options);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
328
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
329 if (nargout == 2)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
330 varargout{1} = t;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
331 varargout{2} = y;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
332 elseif (nargout == 1)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
333 varargout{1}.x = t; # Time stamps are saved in field x
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
334 varargout{1}.y = y; # Results are saved in field y
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
335 varargout{1}.solver = solver;
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
336 if (options.haveeventfunction)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
337 varargout{1}.xe = te; # Time info when an event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
338 varargout{1}.ye = ye; # Results when an event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
339 varargout{1}.ie = ie; # Index info which event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
340 endif
22938
a88ceac2aa53 Fix undefined return argument for more than 2 outputs from ode15i,ode15s.
Rik <rik@octave.org>
parents: 22936
diff changeset
341 elseif (nargout > 2)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
342 varargout = cell (1,5);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
343 varargout{1} = t;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
344 varargout{2} = y;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
345 if (options.haveeventfunction)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
346 varargout{3} = te; # Time info when an event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
347 varargout{4} = ye; # Results when an event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
348 varargout{5} = ie; # Index info which event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
349 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
350 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
351
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
352 endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
353
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
354 function res = wrap (t, y, yp, Mass, havetimedep, havestatedep, fun)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
355
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
356 if (! isempty (Mass) && havestatedep)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
357 res = Mass (t, y) * yp - fun (t, y);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
358 elseif (! isempty (Mass) && havetimedep)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
359 res = Mass (t) * yp - fun (t, y);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
360 elseif (! isempty (Mass))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
361 res = Mass * yp - fun (t, y);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
362 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
363 res = yp - fun (t, y);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
364 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
365
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
366 endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
367
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
368 function [jac, jact] = wrapjacfun (t, y, yp, Jac, Mass,
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
369 havetimedep, havejacfun)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
370 if (havejacfun)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
371 jac = - Jac (t, y);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
372 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
373 jac = - Jac;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
374 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
375
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
376 if (! isempty (Mass) && havetimedep)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
377 jact = Mass (t);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
378 elseif (! isempty (Mass))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
379 jact = Mass;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
380 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
381 jact = speye (numel (y));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
382 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
383
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
384 endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
385
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
386
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
387 %!demo
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
388 %! ## Solve Robertson's equations with ode15s
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
389 %! fun = @ (t, y) [-0.04*y(1) + 1e4*y(2).*y(3);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
390 %! 0.04*y(1) - 1e4*y(2).*y(3) - 3e7*y(2).^2;
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
391 %! y(1) + y(2) + y(3) - 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
392 %!
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
393 %! y0 = [1; 0; 0];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
394 %! tspan = [0, 4*logspace(-6, 6)];
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
395 %! M = [1, 0, 0; 0, 1, 0; 0, 0, 0];
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
396 %!
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
397 %! options = odeset ("RelTol", 1e-4, "AbsTol", [1e-6, 1e-10, 1e-6],
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
398 %! "MStateDependence", "none", "Mass", M);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
399 %!
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
400 %! [t, y] = ode15s (fun, tspan, y0, options);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
401 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
402 %! y(:,2) = 1e4 * y(:,2);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
403 %! figure (2);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
404 %! semilogx (t, y, "o");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
405 %! xlabel ("time");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
406 %! ylabel ("species concentration");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
407 %! title ("Robertson DAE problem with a Conservation Law");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
408 %! legend ("y1", "y2", "y3");
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
409
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
410 %!function ydot = fpol (t, y) # Van der Pol equation
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
411 %! ydot = [y(2); (1 - y(1)^2) * y(2) - y(1)];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
412 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
413 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
414 %!function ref = fref () # The computed reference sol
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
415 %! ref = [0.32331666704577, -1.83297456798624];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
416 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
417 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
418 %!function jac = fjac (t, y) # its Jacobian
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
419 %! jac = [0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
420 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
421 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
422 %!function jac = fjcc (t, y) # sparse type
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
423 %! jac = sparse ([0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
424 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
425 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
426 %!function mas = fmas (t, y)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
427 %! mas = [1, 0; 0, 1]; # Dummy mass matrix for tests
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
428 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
429 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
430 %!function mas = fmsa (t, y)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
431 %! mas = sparse ([1, 0; 0, 1]); # A sparse dummy matrix
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
432 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
433 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
434 %!function res = rob (t, y)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
435 %! res = [-0.04*y(1) + 1e4*y(2).*y(3);
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
436 %! 0.04*y(1) - 1e4*y(2).*y(3) - 3e7*y(2).^2;
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
437 %! y(1) + y(2) + y(3) - 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
438 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
439 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
440 %!function refrob = frefrob ()
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
441 %! refrob = [100, 0.617234887614937, 0.000006153591397, 0.382758958793666];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
442 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
443 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
444 %!function [val, isterminal, direction] = feve (t, y)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
445 %! isterminal = [0, 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
446 %! if (t < 1e1)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
447 %! val = [-1, -2];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
448 %! else
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
449 %! val = [1, 3];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
450 %! endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
451 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
452 %! direction = [1, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
453 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
454 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
455 %!function masrob = massdensefunstate (t, y)
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
456 %! masrob = [1, 0, 0; 0, 1, 0; 0, 0, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
457 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
458 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
459 %!function masrob = masssparsefunstate (t, y)
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
460 %! masrob = sparse ([1, 0, 0; 0, 1, 0; 0, 0, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
461 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
462 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
463 %!function masrob = massdensefuntime (t)
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
464 %! masrob = [1, 0, 0; 0, 1, 0; 0, 0, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
465 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
466 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
467 %!function masrob = masssparsefuntime (t)
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
468 %! masrob = sparse ([1, 0, 0; 0, 1, 0; 0, 0, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
469 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
470 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
471 %!function jac = jacfundense (t, y)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
472 %! jac = [-0.04, 1e4*y(3), 1e4*y(2);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
473 %! 0.04, -1e4*y(3)-6e7*y(2), -1e4*y(2);
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
474 %! 1, 1, 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
475 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
476 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
477 %!function jac = jacfunsparse (t, y)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
478 %! jac = sparse([-0.04, 1e4*y(3), 1e4*y(2);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
479 %! 0.04, -1e4*y(3)-6e7*y(2), -1e4*y(2);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
480 %! 1, 1, 1]);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
481 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
482
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
483 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
484 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
485 %! "Mass", [1, 0, 0; 0, 1, 0; 0, 0, 0]);
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
486 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
487 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
488
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
489 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
490 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
491 %! "Mass", sparse ([1, 0, 0; 0, 1, 0; 0, 0, 0]));
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
492 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
493 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
494
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
495 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
496 %! opt = odeset ("MStateDependence", "none",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
497 %! "Mass", @massdensefunstate);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
498 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
499 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
500
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
501 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
502 %! opt = odeset ("MStateDependence", "none",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
503 %! "Mass", @masssparsefunstate);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
504 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
505 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
506
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
507 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
508 %! opt = odeset ("MStateDependence", "none",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
509 %! "Mass", "massdensefuntime");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
510 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
511 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
512
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
513 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
514 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
515 %! "Mass", [1, 0, 0; 0, 1, 0; 0, 0, 0],
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
516 %! "Jacobian", "jacfundense");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
517 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
518 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
519
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
520 %!testif HAVE_SUNDIALS
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
521 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
522 %! "Mass", sparse ([1, 0, 0; 0, 1, 0; 0, 0, 0]),
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
523 %! "Jacobian", @jacfundense);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
524 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
525 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
526
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
527 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
528 %! warning ("off", "ode15s:mass_state_dependent_provided", "local");
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
529 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
530 %! "Mass", @massdensefunstate,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
531 %! "Jacobian", @jacfundense);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
532 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
533 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
534
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
535 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
536 %! warning ("off", "ode15s:mass_state_dependent_provided", "local");
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
537 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
538 %! "Mass", @masssparsefunstate,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
539 %! "Jacobian", @jacfundense);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
540 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
541 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
542
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
543 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
544 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
545 %! "Mass", @massdensefuntime,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
546 %! "Jacobian", @jacfundense);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
547 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
548 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
549 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
550 %! "Mass", "masssparsefuntime",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
551 %! "Jacobian", "jacfundense");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
552 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
553 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
554
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
555 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
556 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
557 %! "Mass", [1, 0, 0; 0, 1, 0; 0, 0, 0],
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
558 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
559 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
560 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
561
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
562 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
563 %! opt = odeset ("MStateDependence", "none",
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
564 %! "Mass", sparse ([1, 0, 0; 0, 1, 0; 0, 0, 0]),
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
565 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
566 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
567 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
568
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
569 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
570 %! warning ("off", "ode15s:mass_state_dependent_provided", "local");
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
571 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
572 %! "Mass", @massdensefunstate,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
573 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
574 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
575 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
576
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
577 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
578 %! warning ("off", "ode15s:mass_state_dependent_provided", "local");
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
579 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
580 %! "Mass", @masssparsefunstate,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
581 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
582 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
583 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
584
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
585 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
586 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
587 %! "Mass", @massdensefuntime,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
588 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
589 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
590 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
591
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
592 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
593 %! opt = odeset ("MStateDependence", "none",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
594 %! "Mass", @masssparsefuntime,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
595 %! "Jacobian", @jacfunsparse);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
596 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
597 %! assert ([t(end), y(end,:)], frefrob, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
598
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
599 ## two output arguments
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
600 %!testif HAVE_SUNDIALS
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
601 %! [t, y] = ode15s (@fpol, [0, 2], [2, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
602 %! assert ([t(end), y(end,:)], [2, fref], 1e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
603
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
604 ## anonymous function instead of real function
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
605 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
606 %! fvdb = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
607 %! [t, y] = ode15s (fvdb, [0, 2], [2, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
608 %! assert ([t(end), y(end,:)], [2, fref], 1e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
609
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
610 ## Solve another anonymous function below zero
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
611 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
612 %! ref = [0, 14.77810590694212];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
613 %! [t, y] = ode15s (@(t,y) y, [-2, 0], 2);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
614 %! assert ([t(end), y(end,:)], ref, 5e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
615
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
616 ## InitialStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
617 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
618 %! opt = odeset ("InitialStep", 1e-8);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
619 %! [t, y] = ode15s (@fpol, [0, 0.2], [2, 0], opt);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
620 %! assert (t(2)-t(1), 1e-8, 1e-9);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
621
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
622 ## MaxStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
623 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
624 %! opt = odeset ("MaxStep", 1e-3);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
625 %! sol = ode15s (@fpol, [0, 0.2], [2, 0], opt);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
626 %! assert (sol.x(5)-sol.x(4), 1e-3, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
627
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
628 ## Solve in backward direction starting at t=0
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
629 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
630 %! ref = [-1.205364552835178, 0.951542399860817];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
631 %! sol = ode15s (@fpol, [0 -2], [2, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
632 %! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 5e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
633
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
634 ## Solve in backward direction starting at t=2
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
635 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
636 %! ref = [-1.205364552835178, 0.951542399860817];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
637 %! sol = ode15s (@fpol, [2, 0 -2], fref);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
638 %! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 3e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
639
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
640 ## Solve another anonymous function in backward direction
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
641 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
642 %! ref = [-1, 0.367879437558975];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
643 %! sol = ode15s (@(t,y) y, [0 -1], 1);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
644 %! assert ([sol.x(end), sol.y(end,:)], ref, 1e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
645
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
646 ## Solve another anonymous function below zero
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
647 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
648 %! ref = [0, 14.77810590694212];
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
649 %! sol = ode15s (@(t,y) y, [-2, 0], 2);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
650 %! assert ([sol.x(end), sol.y(end,:)], ref, 5e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
651
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
652 ## Solve in backward direction starting at t=0 with MaxStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
653 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
654 %! ref = [-1.205364552835178, 0.951542399860817];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
655 %! opt = odeset ("MaxStep", 1e-3);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
656 %! sol = ode15s (@fpol, [0 -2], [2, 0], opt);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
657 %! assert (abs (sol.x(8)-sol.x(7)), 1e-3, 1e-3);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
658 %! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
659
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
660 ## AbsTol option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
661 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
662 %! opt = odeset ("AbsTol", 1e-5);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
663 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
664 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 4e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
665
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
666 ## AbsTol and RelTol option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
667 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
668 %! opt = odeset ("AbsTol", 1e-8, "RelTol", 1e-8);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
669 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
670 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
671
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
672 ## RelTol option -- higher accuracy
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
673 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
674 %! opt = odeset ("RelTol", 1e-8);
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
675 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
676 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-4);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
677
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
678 ## Mass option as function
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
679 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
680 %! opt = odeset ("Mass", @fmas, "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
681 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
682 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 3e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
683
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
684 ## Mass option as matrix
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
685 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
686 %! opt = odeset ("Mass", eye (2,2), "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
687 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
688 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 3e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
689
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
690 ## Mass option as sparse matrix
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
691 %!testif HAVE_SUNDIALS
23023
63a12df71848 avoid sparse jacobian tests for ode15i and ode15s if IDAKLU is missing
John W. Eaton <jwe@octave.org>
parents: 22941
diff changeset
692 %! opt = odeset ("Mass", speye (2), "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
693 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
694 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 3e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
695
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
696 ## Mass option as function and sparse matrix
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
697 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
698 %! opt = odeset ("Mass", "fmsa", "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
699 %! sol = ode15s (@fpol, [0, 2], [2, 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
700 %! assert ([sol.x(end), sol.y(end,:)], [2, fref], 3e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
701
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
702 ## Refine
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
703 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
704 %! opt2 = odeset ("Refine", 3, "Mass", @massdensefunstate,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
705 %! "MStateDependence", "none");
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
706 %! opt1 = odeset ("Mass", @massdensefunstate, "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
707 %! [t, y] = ode15s (@rob, [0, 100], [1; 0; 0], opt1);
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
708 %! [t2, y2] = ode15s (@rob, [0, 100], [1; 0; 0], opt2);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
709 %! assert (numel (t2), numel (t) * 3, 3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
710
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
711 ## Refine ignored if numel (trange) > 2
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
712 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
713 %! opt2 = odeset ("Refine", 3, "Mass", "massdensefunstate",
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
714 %! "MStateDependence", "none");
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
715 %! opt1 = odeset ("Mass", @massdensefunstate, "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
716 %! [t, y] = ode15s ("rob", [0, 10, 100], [1; 0; 0], opt1);
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
717 %! [t2, y2] = ode15s ("rob", [0, 10, 100], [1; 0; 0], opt2);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22928
diff changeset
718 %! assert (numel (t2), numel (t));
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
719
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
720 ## Events option add further elements in sol
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
721 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
722 %! opt = odeset ("Events", @feve, "Mass", @massdensefunstate,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
723 %! "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
724 %! sol = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
725 %! assert (isfield (sol, "ie"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
726 %! assert (sol.ie, [0;1]);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
727 %! assert (isfield (sol, "xe"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
728 %! assert (isfield (sol, "ye"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
729 %! assert (sol.x(end), 10, 1);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
730
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
731 ## Events option, five output arguments
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
732 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
733 %! opt = odeset ("Events", @feve, "Mass", @massdensefunstate,
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
734 %! "MStateDependence", "none");
22936
f2a1fc90a903 * ode15s.m: More minor style fixes in tests.
John W. Eaton <jwe@octave.org>
parents: 22935
diff changeset
735 %! [t, y, te, ye, ie] = ode15s (@rob, [0, 100], [1; 0; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
736 %! assert ([t(end), te', ie'], [10, 10, 10, 0, 1], [1, 0.5, 0.5, 0, 0]);
24008
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
737
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
738 ## Initial solution as row vector
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
739 %!testif HAVE_SUNDIALS
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
740 %! A = zeros (2);
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
741 %! [tout, yout] = ode15s (@(t, y) A * y, [0, 1], [1, 1]);
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
742 %! assert (yout, ones (18, 2))
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
743
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
744 %!testif HAVE_SUNDIALS
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
745 %! A = zeros (2);
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
746 %! fail ("ode15s (@(t, y) A * y, [0, 1], eye (2))",
24012
3d65720cd68a ode15s.m: Fix message in BIST test (bug #50192).
Rik <rik@octave.org>
parents: 24008
diff changeset
747 %! "ode15s: Y0 must be a numeric vector");