annotate scripts/ode/ode15i.m @ 27923:bd51beb6205e

update formatting of copyright notices * Use <https://octave.org/copyright/> instead of <https://octave.org/COPYRIGHT.html/>. * For consistency with other comments in the Octave sources, use C++-style comments for copyright blocks in C and C++ files. * Use delimiters above and below copyright blocks that are appropriate for the language used in the file. * Eliminate extra spacing inside copyright blocks. * lex.ll (looks_like_copyright): Also allow newlines and carriage returns before the word "Copyright". * scripts/mk-doc.pl (gethelp): Also skip empty comment lines. * bp-table.cc, type.m: Adjust tests.
author John W. Eaton <jwe@octave.org>
date Wed, 08 Jan 2020 11:59:41 -0500
parents 1891570abac8
children c20b7290c778
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: 26908
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: 24008
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: 24008
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: 23086
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: 23086
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: 23086
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: 24008
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}] =} ode15i (@var{fun}, @var{trange}, @var{y0}, @var{yp0})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
28 ## @deftypefnx {} {[@var{t}, @var{y}] =} ode15i (@var{fun}, @var{trange}, @var{y0}, @var{yp0}, @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}] =} ode15i (@dots{})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
30 ## @deftypefnx {} {@var{solution} =} ode15i (@dots{})
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
31 ## @deftypefnx {} {} ode15i (@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 fully-implicit Ordinary Differential Equations (ODEs) or
476fc012d09a doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents: 25054
diff changeset
33 ## 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{ode15i} 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
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
39 ## name of the function that defines the ODE: @code{0 = f(t,y,yp)}. The
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
40 ## function must accept three inputs where the first is time @var{t}, the
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
41 ## second is the function value @var{y} (a column vector), and the third
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
42 ## is the derivative value @var{yp} (a column vector).
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
43 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
44 ## @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
45 ## 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
46 ## 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
47 ## 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
48 ## instances.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
49 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
50 ## @var{y0} and @var{yp0} contain the initial values for the unknowns @var{y}
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
51 ## and @var{yp}. If they are row vectors then the solution @var{y} will be a
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
52 ## matrix in which each column is the solution for the corresponding initial
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
53 ## value in @var{y0} and @var{yp0}.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
54 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
55 ## @var{y0} and @var{yp0} must be consistent initial conditions, meaning that
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
56 ## @code{f(t,y0,yp0) = 0} is satisfied. The function @code{decic} may be used
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
57 ## to compute consistent initial conditions given initial guesses.
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
58 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
59 ## The optional fifth 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
60 ## 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
61 ##
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
62 ## 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
63 ## 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
64 ## 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
65 ## 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
66 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
67 ## 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
68 ## 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
69 ## 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
70 ## 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
71 ## @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
72 ## additional information returned.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
73 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
74 ## 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
75 ## 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
76 ## 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
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 ## 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
79 ## 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
80 ## @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
81 ## 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
82 ## of multiple Event functions.
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
83 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
84 ## Example: Solve @nospell{Robertson's} equations:
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
85 ##
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
86 ## @smallexample
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
87 ## @group
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
88 ## function r = robertson_dae (@var{t}, @var{y}, @var{yp})
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
89 ## r = [ -(@var{yp}(1) + 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
90 ## -(@var{yp}(2) - 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
91 ## @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
92 ## endfunction
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
93 ## [@var{t},@var{y}] = ode15i (@@robertson_dae, [0, 1e3], [1; 0; 0], [-1e-4; 1e-4; 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
94 ## @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
95 ## @end smallexample
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
96 ## @seealso{decic, odeset, odeget}
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
97 ## @end deftypefn
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
98
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
99 function varargout = ode15i (fun, trange, y0, yp0, varargin)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
100
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
101 if (nargin < 4)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
102 print_usage ();
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
103 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
104
25026
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
105 solver = "ode15i";
f886561f9696 doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents: 24534
diff changeset
106
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
107 n = numel (y0);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
108
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
109 if (nargin > 4)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
110 options = varargin{1};
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
111 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
112 options = odeset ();
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
113 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
114
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
115 ## Check fun, trange, y0, yp0
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
116 fun = check_default_input (fun, trange, solver, y0, yp0);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
117
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
118 if (! isempty (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
119 if (ischar (options.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
120 try
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
121 options.Jacobian = str2func (options.Jacobian);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
122 catch
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
123 warning (lasterr);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
124 end_try_catch
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25712
diff changeset
125 if (! is_function_handle (options.Jacobian))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
126 error ("Octave:invalid-input-arg",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
127 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
128 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
129 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
130 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
131
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
132 if (! isempty (options.OutputFcn))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
133 if (ischar (options.OutputFcn))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
134 try
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
135 options.OutputFcn = str2func (options.OutputFcn);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
136 catch
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
137 warning (lasterr);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
138 end_try_catch
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25712
diff changeset
139 if (! is_function_handle (options.OutputFcn))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
140 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
141 [solver ": invalid value assigned to field 'OutputFcn'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
142 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
143 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
144 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
145
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
146 if (! isempty (options.Events))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
147 if (ischar (options.Events))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
148 try
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
149 options.Events = str2func (options.Events);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
150 catch
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
151 warning (lasterr);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
152 end_try_catch
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25712
diff changeset
153 if (! is_function_handle (options.Events)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
154 && ! ismatrix (options.Events))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
155 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
156 [solver ": invalid value assigned to field 'Events'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
157 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
158 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
159 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
160
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
161 [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
162
22902
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
163 persistent ignorefields = {"NonNegative", "Mass", ...
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
164 "MStateDependence", "MvPattern", ...
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
165 "MassSingular", "InitialSlope", "BDF"};
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
166
22902
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
167 defaults = rmfield (defaults, ignorefields);
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
168 classes = rmfield (classes, ignorefields);
284bbb0328f2 Update odei5{i,s} code for Octave 4.3.0+ and Sundials 2.7.0.
Carlo de Falco <carlo.defalco@polimi.it>
parents: 22901
diff changeset
169 attributes = rmfield (attributes, ignorefields);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
170
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
171 classes = odeset (classes, "Vectorized", {});
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
172 attributes = odeset (attributes, "Jacobian", {}, "Vectorized", {});
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
173
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
174 options = odemergeopts ("ode15i", options, defaults,
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
175 classes, attributes, solver);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
176
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
177 ## Jacobian
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
178 options.havejac = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
179 options.havejacsparse = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
180 options.havejacfun = 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.Jacobian))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
183 options.havejac = true;
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
184 if (iscell (options.Jacobian))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
185 if (numel (options.Jacobian) == 2)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
186 if (issparse (options.Jacobian{1}) && issparse (options.Jacobian{2}))
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
187 options.havejacsparse = true; # Jac is sparse cell
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
188 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
189
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
190 if (any (size (options.Jacobian{1}) != [n n])
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
191 || any (size (options.Jacobian{2}) != [n n])
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
192 || ! isnumeric (options.Jacobian{1})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
193 || ! isnumeric (options.Jacobian{2})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
194 || ! isreal (options.Jacobian{1})
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
195 || ! isreal (options.Jacobian{2}))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
196 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
197 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
198 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
199 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
200 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
201 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
202 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
203
25803
23483673ba43 Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents: 25712
diff changeset
204 elseif (is_function_handle (options.Jacobian))
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
205 options.havejacfun = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
206 if (nargin (options.Jacobian) == 3)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
207 [A, B] = options.Jacobian (trange(1), y0, yp0);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
208 if (issparse (A) && issparse (B))
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
209 options.havejacsparse = true; # Jac is sparse fun
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
210 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
211
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
212 if (any (size (A) != [n n]) || any (size (B) != [n n])
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
213 || ! isnumeric (A) || ! isnumeric (B) || ! isreal (A)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
214 || ! isreal (B))
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
215 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
216 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
217 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
218 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
219 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
220 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
221 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
222 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
223 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
224 [solver ": invalid value assigned to field 'Jacobian'"]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
225 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
226 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
227
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
228 ## Abstol and Reltol
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
229 options.haveabstolvec = false;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
230
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
231 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
232 error ("Octave:invalid-input-arg",
22941
d92ec2901770 Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents: 22938
diff changeset
233 [solver ": invalid value assigned to field 'AbsTol'"]);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
234
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
235 elseif (numel (options.AbsTol) == n)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
236 options.haveabstolvec = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
237 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
238
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
239 ## Stats
22935
c9344df03da5 Allow case-insensitive option argument 'on' to ode solvers (bug #49918).
Rik <rik@octave.org>
parents: 22933
diff changeset
240 options.havestats = strcmpi (options.Stats, "on");
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
241
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
242 ## 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
243 if (nargout == 1)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
244 options.Refine = 1;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
245 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
246
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
247 ## OutputFcn and OutputSel
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
248 if (isempty (options.OutputFcn) && nargout == 0)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
249 options.OutputFcn = @odeplot;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
250 options.haveoutputfunction = true;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
251 else
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
252 options.haveoutputfunction = ! isempty (options.OutputFcn);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
253 endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
254
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
255 options.haveoutputselection = ! isempty (options.OutputSel);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
256 if (options.haveoutputselection)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
257 options.OutputSel = options.OutputSel - 1;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
258 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
259
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
260 ## Events
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
261 options.haveeventfunction = ! isempty (options.Events);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
262
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
263 [t, y, te, ye, ie] = __ode15__ (fun, trange, y0, yp0, options);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
264
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
265 if (nargout == 2)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
266 varargout{1} = t;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
267 varargout{2} = y;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
268 elseif (nargout == 1)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
269 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: 22931
diff changeset
270 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
271 varargout{1}.solver = solver;
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
272 if (options.haveeventfunction)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
273 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
274 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
275 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
276 endif
22938
a88ceac2aa53 Fix undefined return argument for more than 2 outputs from ode15i,ode15s.
Rik <rik@octave.org>
parents: 22935
diff changeset
277 elseif (nargout > 2)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
278 varargout = cell (1,5);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
279 varargout{1} = t;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
280 varargout{2} = y;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
281 if (options.haveeventfunction)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
282 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
283 varargout{4} = ye; # Results when an event occurred
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
284 varargout{5} = ie; # Index info which event occurred
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 endif
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
287
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
288 endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
289
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
290
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
291 %!demo
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
292 %! ## Solve Robertson's equations with ode15i
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
293 %! fun = @(t, y, yp) [-(yp(1) + 0.04*y(1) - 1e4*y(2)*y(3));
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
294 %! -(yp(2) - 0.04*y(1) + 1e4*y(2)*y(3) + 3e7*y(2)^2);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
295 %! y(1) + y(2) + y(3) - 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
296 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
297 %! opt = odeset ("RelTol", 1e-4, "AbsTol", [1e-8, 1e-14, 1e-6]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
298 %! y0 = [1; 0; 0];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
299 %! yp0 = [-1e-4; 1e-4; 0];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
300 %! tspan = [0 4*logspace(-6, 6)];
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
301 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
302 %! [t, y] = ode15i (fun, tspan, y0, yp0, opt);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
303 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
304 %! y(:,2) = 1e4 * y(:, 2);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
305 %! figure (2);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
306 %! semilogx (t, y, "o");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
307 %! xlabel ("time");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
308 %! ylabel ("species concentration");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
309 %! title ("Robertson DAE problem with a Conservation Law");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
310 %! legend ("y1", "y2", "y3");
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
311
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
312 %!function res = rob (t, y, yp)
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
313 %! res =[-(yp(1) + 0.04*y(1) - 1e4*y(2)*y(3));
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
314 %! -(yp(2) - 0.04*y(1) + 1e4*y(2)*y(3) + 3e7*y(2)^2);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
315 %! y(1) + y(2) + y(3) - 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
316 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
317 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
318 %!function ref = fref ()
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
319 %! ref = [100, 0.617234887614937, 0.000006153591397, 0.382758958793666];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
320 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
321 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
322 %!function ref2 = fref2 ()
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
323 %! ref2 = [4e6 0 0 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
324 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
325 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
326 %!function [DFDY, DFDYP] = jacfundense (t, y, yp)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
327 %! DFDY = [-0.04, 1e4*y(3), 1e4*y(2);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
328 %! 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
329 %! 1, 1, 1];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
330 %! DFDYP = [-1, 0, 0;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
331 %! 0, -1, 0;
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
332 %! 0, 0, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
333 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
334 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
335 %!function [DFDY, DFDYP] = jacfunsparse (t, y, yp)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
336 %! DFDY = 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
337 %! 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
338 %! 1, 1, 1]);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
339 %! DFDYP = sparse ([-1, 0, 0;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
340 %! 0, -1, 0;
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
341 %! 0, 0, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
342 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
343 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
344 %!function [DFDY, DFDYP] = jacwrong (t, y, yp)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
345 %! DFDY = [-0.04, 1e4*y(3);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
346 %! 0.04, -1e4*y(3)-6e7*y(2)];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
347 %! DFDYP = [-1, 0;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
348 %! 0, -1];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
349 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
350 %!
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
351 %!function [DFDY, DFDYP, A] = jacwrong2 (t, y, yp)
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
352 %! DFDY = [-0.04, 1e4*y(3), 1e4*y(2);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
353 %! 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
354 %! 1, 1, 1];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
355 %! DFDYP = [-1, 0, 0;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
356 %! 0, -1, 0;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
357 %! 0, 0, 0];
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
358 %! A = DFDY;
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
359 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
360 %!
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
361 %!function [val, isterminal, direction] = ff (t, y, yp)
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
362 %! isterminal = [0, 1];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
363 %! if (t < 1e1)
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
364 %! val = [-1, -2];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
365 %! else
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
366 %! val = [1, 3];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
367 %! endif
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
368 %!
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
369 %! direction = [1, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
370 %!endfunction
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
371
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
372 ## anonymous function instead of real function
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
373 %!testif HAVE_SUNDIALS
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
374 %! ref = 0.049787079136413;
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
375 %! ff = @(t, u, udot) udot + 3 * u;
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
376 %! [t, y] = ode15i (ff, 0:1, 1, -3);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
377 %! assert ([t(end), y(end)], [1, ref], 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
378
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
379 ## function passed as string
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
380 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
381 %! [t, y] = ode15i ("rob", [0, 100, 200], [1; 0; 0], [-1e-4; 1e-4; 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
382 %! assert ([t(2), y(2,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
383
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
384 ## solve in intermidiate step
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
385 %!testif HAVE_SUNDIALS
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
386 %! [t, y] = ode15i (@rob, [0, 100, 200], [1; 0; 0], [-1e-4; 1e-4; 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
387 %! assert ([t(2), y(2,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
388
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
389 ## numel(trange) = 2 final value
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
390 %!testif HAVE_SUNDIALS
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
391 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
392 %! assert ([t(end), y(end,:)], fref, 1e-5);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
393
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
394 ## With empty options
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
395 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
396 %! opt = odeset();
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
397 %! [t, y] = ode15i (@rob, [0, 1e6, 2e6, 3e6, 4e6], [1; 0; 0],
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
398 %! [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
399 %! assert ([t(end), y(end,:)], fref2, 1e-3);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
400 %! opt = odeset();
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
401
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
402 ## Without options
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
403 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
404 %! [t, y] = ode15i (@rob, [0, 1e6, 2e6, 3e6, 4e6], [1; 0; 0],[-1e-4; 1e-4; 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
405 %! assert ([t(end), y(end,:)], fref2, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
406
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
407 ## InitialStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
408 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
409 %! opt = odeset ("InitialStep", 1e-8);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
410 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
411 %! 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
412
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
413 ## MaxStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
414 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
415 %! opt = odeset ("MaxStep", 1e-3);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
416 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
417 %! assert (t(5)-t(4), 1e-3, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
418
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
419 ## AbsTol scalar option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
420 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
421 %! opt = odeset ("AbsTol", 1e-8);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
422 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
423 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
424
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
425 ## AbsTol scalar and RelTol option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
426 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
427 %! opt = odeset ("AbsTol", 1e-8, "RelTol", 1e-6);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
428 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
429 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
430
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
431 ## AbsTol vector option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
432 %!testif HAVE_SUNDIALS
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
433 %! opt = odeset ("AbsTol", [1e-8, 1e-14, 1e-6]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
434 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
435 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
436
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
437 ## AbsTol vector and RelTol option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
438 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
439 %! opt = odeset ("AbsTol", [1e-8, 1e-14,1e-6], "RelTol", 1e-6);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
440 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4;1e-4;0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
441 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
442
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
443 ## RelTol option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
444 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
445 %! opt = odeset ("RelTol", 1e-6);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
446 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
447 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
448
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
449 ## Jacobian fun dense
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
450 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
451 %! opt = odeset ("Jacobian", @jacfundense);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
452 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
453 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
454
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
455 ## Jacobian fun dense as string
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
456 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
457 %! opt = odeset ("Jacobian", "jacfundense");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
458 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
459 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
460
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
461 ## Jacobian fun sparse
26894
ee6300e77c92 Update detection of sundials in the build system (bug #52475).
Carlo de Falco <carlo.defalco@polimi.it>
parents: 26376
diff changeset
462 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
463 %! opt = odeset ("Jacobian", @jacfunsparse, "AbsTol", 1e-7, "RelTol", 1e-7);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
464 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
465 %! assert ([t(end), y(end,:)], fref, 1e-3);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
466
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
467 ## Solve in backward direction starting at t=100
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
468 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
469 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
470 %! Yref = [0.617234887614937, 0.000006153591397, 0.382758958793666];
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
471 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
472 %! [t2, y2] = ode15i (@rob, [100, 0], Yref', YPref);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
473 %! assert ([t2(end), y2(end,:)], [0, 1, 0, 0], 2e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
474
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
475 ## Solve in backward direction with MaxStep option
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
476 #%!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
477 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
478 %! Yref = [0.617234887614937, 0.000006153591397, 0.382758958793666];
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
479 %! opt = odeset ("MaxStep", 1e-2);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
480 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
481 %! [t2, y2] = ode15i (@rob, [100, 0], Yref', YPref, opt);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
482 %! assert ([t2(end), y2(end,:)], [0, 1, 0, 0], 2e-2);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
483 %! assert (t2(9)-t2(10), 1e-2, 1e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
484
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
485 ## Solve in backward direction starting with intermediate step
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
486 #%!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
487 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
488 %! Yref = [0.617234887614937, 0.000006153591397, 0.382758958793666];
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
489 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
490 %! [t2, y2] = ode15i (@rob, [100, 5, 0], Yref', YPref);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
491 %! assert ([t2(end), y2(end,:)], [0, 1, 0, 0], 2e-2);
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
492
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
493 ## Refine
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
494 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
495 %! opt = odeset ("Refine", 3);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
496 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
497 %! [t2, y2] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
498 %! 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
499
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
500 ## Refine ignored if numel (trange) > 2
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
501 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
502 %! opt = odeset ("Refine", 3);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
503 %! [t, y] = ode15i (@rob, [0, 10, 100], [1; 0; 0], [-1e-4; 1e-4; 0]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
504 %! [t2, y2] = ode15i (@rob, [0, 10, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
505 %! assert (numel (t2), numel (t));
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 ## 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
508 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
509 %! opt = odeset ("Events", @ff);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
510 %! sol = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0], opt);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
511 %! assert (isfield (sol, "ie"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
512 %! assert (sol.ie, [0;1]);
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
513 %! assert (isfield (sol, "xe"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
514 %! assert (isfield (sol, "ye"));
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
515 %! 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
516
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
517 ## Events option, five output arguments
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
518 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
519 %! opt = odeset ("Events", @ff);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
520 %! [t, y, te, ye, ie] = ode15i (@rob, [0, 100], [1; 0; 0],
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
521 %! [-1e-4; 1e-4; 0], opt);
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
522 %! assert ([t(end), te', ie'], [10, 10, 10, 0, 1], [1, 0.2, 0.2, 0, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
523
24008
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
524 ## Initial solutions as row vectors
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
525 %!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
526 %! A = eye (2);
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
527 %! [tout, yout] = ode15i (@(t, y, yp) A * y - A * yp, ...
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
528 %! [0, 1], [1, 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
529 %! assert (size (yout), [20, 2])
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
530
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
531 %!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
532 %! A = eye (2);
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
533 %! [tout, yout] = ode15i (@(t, y, yp) A * y - A * yp, ...
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
534 %! [0, 1], [1, 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
535 %! assert (size (yout), [20, 2])
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
536
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
537 ## Jacobian fun wrong dimension
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
538 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
539 %! opt = odeset ("Jacobian", @jacwrong);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
540 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
541 %! "ode15i: invalid value assigned to field 'Jacobian'");
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 ## Jacobian cell dense wrong dimension
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
544 %!testif HAVE_SUNDIALS
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
545 %! DFDY = [-0.04, 1;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
546 %! 0.04, 1];
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
547 %! DFDYP = [-1, 0, 0;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
548 %! 0, -1, 0;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
549 %! 0, 0, 0];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
550 %! opt = odeset ("Jacobian", {DFDY, DFDYP});
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
551 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
552 %! "invalid value assigned to field 'Jacobian'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
553
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
554 ## Jacobian cell sparse wrong dimension
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
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
556 %! DFDY = sparse ([-0.04, 1;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
557 %! 0.04, 1]);
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
558 %! DFDYP = sparse ([-1, 0, 0;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
559 %! 0, -1, 0;
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
560 %! 0, 0, 0]);
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
561 %! opt = odeset ("Jacobian", {DFDY, DFDYP});
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
562 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
563 %! "invalid value assigned to field 'Jacobian'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
564
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
565 ## Jacobian cell wrong number of matrices
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
566 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
567 %! A = [1 2 3; 4 5 6; 7 8 9];
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
568 %! opt = odeset ("Jacobian", {A,A,A});
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
569 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
570 %! "invalid value assigned to field 'Jacobian'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
571
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
572 ## Jacobian single matrix
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
573 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
574 %! opt = odeset ("Jacobian", [1 2 3; 4 5 6; 7 8 9]);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
575 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
576 %! "invalid value assigned to field 'Jacobian'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
577
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
578 ## Jacobian single matrix wrong dimension
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
579 %!testif HAVE_SUNDIALS
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
580 %! opt = odeset ("Jacobian", [1 2 3; 4 5 6]);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
581 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
582 %! "invalid value assigned to field 'Jacobian'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
583
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
584 ## Jacobian strange field
26908
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
585 ## FIXME: for compatibility with Matlab, it is no longer an error to
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
586 ## create a handle to a nonexistent function and ode15i just checks that
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
587 ## the argument is a function handle (it is) so there is no error for
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
588 ## this case now. If this check really must be done early, before
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
589 ## attempting to use the invalid function handle, then one way is to
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
590 ## call "functions" on it and catch the error if it is invalid.
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
591 ##%!testif HAVE_SUNDIALS
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
592 ##%! opt = odeset ("Jacobian", "_5yVNhWVJWJn47RKnzxPsyb_");
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
593 ##%! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)",
6bd408bde569 * ode15i.m: Comment out failing test
John W. Eaton <jwe@octave.org>
parents: 26894
diff changeset
594 ##%! "invalid value assigned to field 'Jacobian'");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
595
22910
23847979b91e maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents: 22902
diff changeset
596 %!function ydot = fun (t, y, yp)
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
597 %! ydot = [y - yp];
22901
4c56f3ffec04 Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff changeset
598 %!endfunction
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
599
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
600 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
601 %! fail ("ode15i ()", "Invalid call to ode15i");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
602
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
603 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
604 %! fail ("ode15i (1)", "Invalid call to ode15i");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
605
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
606 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
607 %! fail ("ode15i (1, 1)", "Invalid call to ode15i");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
608
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
609 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
610 %! fail ("ode15i (1, 1, 1)", "Invalid call to ode15i");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
611
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
612 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
613 %! fail ("ode15i (1, 1, 1, 1)", "ode15i: fun must be of class:");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
614
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
615 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
616 %! fail ("ode15i (1, 1, 1, 1, 1)", "ode15i: fun must be of class:");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
617
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
618 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
619 %! fail ("ode15i (1, 1, 1, 1, 1, 1)", "ode15i: fun must be of class:");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
620
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
621 %!testif HAVE_SUNDIALS
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
622 %! fail ("ode15i (@fun, 1, 1, 1)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
623 %! "ode15i: invalid value assigned to field 'trange'");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
624
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
625 %!testif HAVE_SUNDIALS
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
626 %! fail ("ode15i (@fun, [1, 1], 1, 1)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
627 %! "ode15i: invalid value assigned to field 'trange'");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
628
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
629 %!testif HAVE_SUNDIALS
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
630 %! fail ("ode15i (@fun, [1, 2], 1, [1, 2])",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
631 %! "ode15i: y0 must have 2 elements");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
632
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
633 %!testif HAVE_SUNDIALS
25712
5625b2237a4d * ode15i.m: In tests, use unlikely symbol name instead of "foo".
John W. Eaton <jwe@octave.org>
parents: 25113
diff changeset
634 %! opt = odeset ("RelTol", "_5yVNhWVJWJn47RKnzxPsyb_");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
635 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
636 %! "ode15i: RelTol must be of class:");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
637
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
638 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
639 %! opt = odeset ("RelTol", [1, 2]);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
640 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
641 %! "ode15i: RelTol must be scalar");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
642
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
643 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
644 %! opt = odeset ("RelTol", -2);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
645 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
646 %! "ode15i: RelTol must be positive");
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
647
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
648 %!testif HAVE_SUNDIALS
25712
5625b2237a4d * ode15i.m: In tests, use unlikely symbol name instead of "foo".
John W. Eaton <jwe@octave.org>
parents: 25113
diff changeset
649 %! opt = odeset ("AbsTol", "_5yVNhWVJWJn47RKnzxPsyb_");
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
650 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
651 %! "ode15i: AbsTol must be of class:");
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
652
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
653 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
654 %! opt = odeset ("AbsTol", -1);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
655 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
656 %! "ode15i: AbsTol must be positive");
22928
dec22bceafa2 use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents: 22910
diff changeset
657
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
658 %!testif HAVE_SUNDIALS
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
659 %! opt = odeset ("AbsTol", [1, 1, 1]);
22931
8133da976602 * ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents: 22928
diff changeset
660 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)",
22933
c3428bb9aca9 ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22931
diff changeset
661 %! "ode15i: invalid value assigned to field 'AbsTol'");
24008
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
662
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
663 %!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
664 %! 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
665 %! fail ("ode15i (@(t, y, yp) A * y - A * yp, [0, 1], eye (2), [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
666 %! "ode15i: Y0 must be a numeric vector");
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
667
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
668 %!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
669 %! 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
670 %! fail ("ode15i (@(t, y, yp) A * y - A * yp, [0, 1], [1, 1], eye (2))",
6e7bb85e32b8 Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents: 24007
diff changeset
671 %! "ode15i: YP0 must be a numeric vector");