Mercurial > octave
annotate scripts/ode/ode15i.m @ 29296:7bf91e98bfc6
__ode15__(): Consider the correct number of input arguments in the event callback (bug #59477).
* ode15i.m: Provide three input arguments (t, y, yp) in the event callback.
* ode15s.m: Provide two input arguments (t, y) in the event callback.
* __ode15__.cc: Provide the number of input arguments in the event callback.
author | Markus Meisinger <chloros2@gmx.de> |
---|---|
date | Mon, 11 Jan 2021 20:11:45 +0100 |
parents | d25cd81b0577 |
children | 7854d5752dd2 |
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 ## |
28713
28d2511f2af2
maint: grammarcheck documentation ahead of 6.1 release.
Rik <rik@octave.org>
parents:
28660
diff
changeset
|
74 ## If no output arguments are requested, and no @qcode{"OutputFcn"} is |
28d2511f2af2
maint: grammarcheck documentation ahead of 6.1 release.
Rik <rik@octave.org>
parents:
28660
diff
changeset
|
75 ## specified in @var{ode_opt}, then the @qcode{"OutputFcn"} is set to |
28d2511f2af2
maint: grammarcheck documentation ahead of 6.1 release.
Rik <rik@octave.org>
parents:
28660
diff
changeset
|
76 ## @code{odeplot} and the 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)) |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
120 if (! exist (options.Jacobian)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
121 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
122 ['ode15i: "Jacobian" function "' options.Jacobian '" not found']); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
123 endif |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
124 options.Jacobian = str2func (options.Jacobian); |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
125 endif |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
126 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
127 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
128 if (! isempty (options.OutputFcn)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
129 if (ischar (options.OutputFcn)) |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
130 if (! exist (options.OutputFcn)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
131 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
132 ['ode15i: "OutputFcn" function "' options.OutputFcn '" not found']); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
133 endif |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
134 options.OutputFcn = str2func (options.OutputFcn); |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
135 endif |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
136 if (! is_function_handle (options.OutputFcn)) |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
137 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
138 'ode15i: "OutputFcn" must be a valid function handle'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
139 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
140 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
141 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
142 if (! isempty (options.Events)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
143 if (ischar (options.Events)) |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
144 if (! exist (options.Events)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
145 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
146 ['ode15i: "Events" function "' options.Events '" not found']); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
147 endif |
28041
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
148 options.Events = str2func (options.Events); |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
149 endif |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
150 if (! is_function_handle (options.Events)) |
5e44268dca6f
Replace input validation relying on str2func with alternatives (bug #57351).
Rik <rik@octave.org>
parents:
27957
diff
changeset
|
151 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
152 'ode15i: "Events" must be a valid function handle'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
153 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
154 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
155 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
156 [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
|
157 |
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
|
158 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
|
159 "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
|
160 "MassSingular", "InitialSlope", "BDF"}; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
161 |
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
|
162 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
|
163 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
|
164 attributes = rmfield (attributes, ignorefields); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
165 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
166 classes = odeset (classes, "Vectorized", {}); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
167 attributes = odeset (attributes, "Jacobian", {}, "Vectorized", {}); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
168 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
169 options = odemergeopts ("ode15i", options, defaults, |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
170 classes, attributes, solver); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
171 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
172 ## Jacobian |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
173 options.havejac = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
174 options.havejacsparse = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
175 options.havejacfun = false; |
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 if (! isempty (options.Jacobian)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
178 options.havejac = true; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
179 if (iscell (options.Jacobian)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
180 if (numel (options.Jacobian) == 2) |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
181 J1 = options.Jacobian{1}; |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
182 J2 = options.Jacobian{2}; |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
183 if ( ! issquare (J1) || ! issquare (J2) |
28660
dc80e087df4b
maint: Strip trailing spaces from files.
Rik <rik@octave.org>
parents:
28047
diff
changeset
|
184 || rows (J1) != n || rows (J2) != n |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
185 || ! isnumeric (J1) || ! isnumeric (J2) |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
186 || ! isreal (J1) || ! isreal (J2)) |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
187 error ("Octave:invalid-input-arg", |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
188 'ode15i: "Jacobian" matrices must be real square matrices'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
189 endif |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
190 if (issparse (J1) && issparse (J2)) |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
191 options.havejacsparse = true; # Jac is sparse cell |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
192 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
193 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
194 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
195 'ode15i: invalid value assigned to field "Jacobian"'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
196 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
197 |
25803
23483673ba43
Use is_function_handle instead of isa (x, "function_handle").
Rik <rik@octave.org>
parents:
25712
diff
changeset
|
198 elseif (is_function_handle (options.Jacobian)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
199 options.havejacfun = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
200 if (nargin (options.Jacobian) == 3) |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
201 [J1, J2] = options.Jacobian (trange(1), y0, yp0); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
202 |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
203 if ( ! issquare (J1) || rows (J1) != n |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
204 || ! isnumeric (J1) || ! isreal (J1) |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
205 || ! issquare (J2) || rows (J2) != n |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
206 || ! isnumeric (J2) || ! isreal (J2)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
207 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
208 'ode15i: "Jacobian" function must evaluate to a real square matrix'); |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
209 endif |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
210 if (issparse (J1) && issparse (J2)) |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
211 options.havejacsparse = true; # Jac is sparse fun |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
212 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
213 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
214 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
215 'ode15i: invalid value assigned to field "Jacobian"'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
216 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
217 else |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
218 error ("Octave:invalid-input-arg", |
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
219 'ode15i: "Jacobian" field must be a function handle or 2-element cell array of square matrices'); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
220 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
221 endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
222 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
223 ## Abstol and Reltol |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
224 options.haveabstolvec = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
225 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
226 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
|
227 error ("Octave:invalid-input-arg", |
28045
13dba3c069f8
Update input validation for odeXXX.m functions.
Rik <rik@octave.org>
parents:
28044
diff
changeset
|
228 'ode15i: 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
|
229 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
230 elseif (numel (options.AbsTol) == n) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
231 options.haveabstolvec = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
232 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
233 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
234 ## Stats |
22935
c9344df03da5
Allow case-insensitive option argument 'on' to ode solvers (bug #49918).
Rik <rik@octave.org>
parents:
22933
diff
changeset
|
235 options.havestats = strcmpi (options.Stats, "on"); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
236 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
237 ## 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
|
238 if (nargout == 1) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
239 options.Refine = 1; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
240 endif |
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 ## OutputFcn and OutputSel |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
243 if (isempty (options.OutputFcn) && nargout == 0) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
244 options.OutputFcn = @odeplot; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
245 options.haveoutputfunction = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
246 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
247 options.haveoutputfunction = ! isempty (options.OutputFcn); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
248 endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
249 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
250 options.haveoutputselection = ! isempty (options.OutputSel); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
251 if (options.haveoutputselection) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
252 options.OutputSel = options.OutputSel - 1; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
253 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
254 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
255 ## Events |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
256 options.haveeventfunction = ! isempty (options.Events); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
257 |
29296
7bf91e98bfc6
__ode15__(): Consider the correct number of input arguments in the event callback (bug #59477).
Markus Meisinger <chloros2@gmx.de>
parents:
29042
diff
changeset
|
258 ## 3 arguments in the event callback of ode15i |
7bf91e98bfc6
__ode15__(): Consider the correct number of input arguments in the event callback (bug #59477).
Markus Meisinger <chloros2@gmx.de>
parents:
29042
diff
changeset
|
259 [t, y, te, ye, ie] = __ode15__ (fun, trange, y0, yp0, options, 3); |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
260 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
261 if (nargout == 2) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
262 varargout{1} = t; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
263 varargout{2} = y; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
264 elseif (nargout == 1) |
29025
7f103819617d
ode15i.m, ode15s.m: transpose struct output for Matlab compatibility (bug #59416)
Rik <rik@octave.org>
parents:
28929
diff
changeset
|
265 varargout{1}.x = t.'; # Time stamps saved in field x (row vector) |
7f103819617d
ode15i.m, ode15s.m: transpose struct output for Matlab compatibility (bug #59416)
Rik <rik@octave.org>
parents:
28929
diff
changeset
|
266 varargout{1}.y = y.'; # Results are saved in field y (row vector) |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
267 varargout{1}.solver = solver; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
268 if (options.haveeventfunction) |
29041
2190720bca3e
ode15i.m, ode15s.m: Transpose event outputs when single struct output requested (bug #59416).
Rik <rik@octave.org>
parents:
29039
diff
changeset
|
269 varargout{1}.xe = te.'; # Time info when an event occurred |
2190720bca3e
ode15i.m, ode15s.m: Transpose event outputs when single struct output requested (bug #59416).
Rik <rik@octave.org>
parents:
29039
diff
changeset
|
270 varargout{1}.ye = ye.'; # Results when an event occurred |
2190720bca3e
ode15i.m, ode15s.m: Transpose event outputs when single struct output requested (bug #59416).
Rik <rik@octave.org>
parents:
29039
diff
changeset
|
271 varargout{1}.ie = ie.'; # Index info which event occurred |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
272 endif |
22938
a88ceac2aa53
Fix undefined return argument for more than 2 outputs from ode15i,ode15s.
Rik <rik@octave.org>
parents:
22935
diff
changeset
|
273 elseif (nargout > 2) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
274 varargout = cell (1,5); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
275 varargout{1} = t; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
276 varargout{2} = y; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
277 if (options.haveeventfunction) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
278 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
|
279 varargout{4} = ye; # Results when an event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
280 varargout{5} = ie; # Index info which event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
281 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
282 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
283 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
284 endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
285 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
286 |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
287 %!demo |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
288 %! ## 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
|
289 %! 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
|
290 %! -(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
|
291 %! y(1) + y(2) + y(3) - 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
292 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
293 %! 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
|
294 %! y0 = [1; 0; 0]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
295 %! yp0 = [-1e-4; 1e-4; 0]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
296 %! 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
|
297 %! |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
298 %! [t, y] = ode15i (fun, tspan, y0, yp0, opt); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
299 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
300 %! y(:,2) = 1e4 * y(:, 2); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
301 %! figure (2); |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
302 %! semilogx (t, y, "o"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
303 %! xlabel ("time"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
304 %! ylabel ("species concentration"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
305 %! 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
|
306 %! legend ("y1", "y2", "y3"); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
307 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
308 %!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
|
309 %! 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
|
310 %! -(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
|
311 %! y(1) + y(2) + y(3) - 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
312 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
313 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
314 %!function ref = fref () |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
315 %! 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
|
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 ref2 = fref2 () |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
319 %! ref2 = [4e6 0 0 1]; |
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 [DFDY, DFDYP] = jacfundense (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
323 %! 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
|
324 %! 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
|
325 %! 1, 1, 1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
326 %! DFDYP = [-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
327 %! 0, -1, 0; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
328 %! 0, 0, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
329 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
330 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
331 %!function [DFDY, DFDYP] = jacfunsparse (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
332 %! 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
|
333 %! 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
|
334 %! 1, 1, 1]); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
335 %! DFDYP = sparse ([-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
336 %! 0, -1, 0; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
337 %! 0, 0, 0]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
338 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
339 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
340 %!function [DFDY, DFDYP] = jacwrong (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
341 %! DFDY = [-0.04, 1e4*y(3); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
342 %! 0.04, -1e4*y(3)-6e7*y(2)]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
343 %! DFDYP = [-1, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
344 %! 0, -1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
345 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
346 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
347 %!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
|
348 %! 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
|
349 %! 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
|
350 %! 1, 1, 1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
351 %! DFDYP = [-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
352 %! 0, -1, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
353 %! 0, 0, 0]; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
354 %! A = DFDY; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
355 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
356 %! |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
357 %!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
|
358 %! isterminal = [0, 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
359 %! if (t < 1e1) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
360 %! val = [-1, -2]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
361 %! else |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
362 %! val = [1, 3]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
363 %! endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
364 %! |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
365 %! direction = [1, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
366 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
367 |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
368 ## anonymous function instead of real function |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
369 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
370 %! ref = 0.049787079136413; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
371 %! ff = @(t, u, udot) udot + 3 * u; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
372 %! [t, y] = ode15i (ff, 0:1, 1, -3); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
373 %! 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
|
374 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
375 ## function passed as string |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
376 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
377 %! [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
|
378 %! 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
|
379 |
27956
2310164737b3
fix many spelling errors (bug #57613)
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
380 ## solve in intermediate step |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
381 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
382 %! [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
|
383 %! 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
|
384 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
385 ## numel(trange) = 2 final value |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
386 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
387 %! [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
|
388 %! 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
|
389 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
390 ## With empty options |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
391 %!testif HAVE_SUNDIALS |
28912
0de38a6ef693
maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents:
28713
diff
changeset
|
392 %! opt = odeset (); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
393 %! [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
|
394 %! [-1e-4; 1e-4; 0], opt); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
395 %! assert ([t(end), y(end,:)], fref2, 1e-3); |
28912
0de38a6ef693
maint: Use Octave convention of space after function name in scripts dir.
Rik <rik@octave.org>
parents:
28713
diff
changeset
|
396 %! opt = odeset (); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
397 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
398 ## Without options |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
399 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
400 %! [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
|
401 %! 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
|
402 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
403 ## InitialStep option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
404 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
405 %! 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
|
406 %! [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
|
407 %! 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
|
408 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
409 ## MaxStep option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
410 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
411 %! 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
|
412 %! [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
|
413 %! 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
|
414 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
415 ## AbsTol scalar option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
416 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
417 %! 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
|
418 %! [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
|
419 %! 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
|
420 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
421 ## AbsTol scalar and RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
422 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
423 %! 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
|
424 %! [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
|
425 %! 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
|
426 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
427 ## AbsTol vector option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
428 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
429 %! 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
|
430 %! [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
|
431 %! 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
|
432 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
433 ## AbsTol vector and RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
434 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
435 %! 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
|
436 %! [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
|
437 %! 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
|
438 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
439 ## RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
440 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
441 %! 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
|
442 %! [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
|
443 %! 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
|
444 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
445 ## Jacobian fun dense |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
446 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
447 %! 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
|
448 %! [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
|
449 %! 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
|
450 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
451 ## Jacobian fun dense as string |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
452 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
453 %! 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
|
454 %! [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
|
455 %! 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
|
456 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
457 ## 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
|
458 %!testif HAVE_SUNDIALS_SUNLINSOL_KLU |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
459 %! 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
|
460 %! [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
|
461 %! 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
|
462 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
463 ## 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
|
464 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
465 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
466 %! 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
|
467 %! [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
|
468 %! [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
|
469 %! 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
|
470 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
471 ## 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
|
472 #%!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
473 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
474 %! 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
|
475 %! 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
|
476 %! [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
|
477 %! [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
|
478 %! 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
|
479 %! 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
|
480 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
481 ## 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
|
482 #%!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
483 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
484 %! 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
|
485 %! [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
|
486 %! [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
|
487 %! 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
|
488 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
489 ## Refine |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
490 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
491 %! 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
|
492 %! [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
|
493 %! [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
|
494 %! 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
|
495 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
496 ## Refine ignored if numel (trange) > 2 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
497 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
498 %! 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
|
499 %! [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
|
500 %! [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
|
501 %! assert (numel (t2), numel (t)); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
502 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
503 ## 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
|
504 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
505 %! 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
|
506 %! 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
|
507 %! assert (isfield (sol, "ie")); |
29042
d25cd81b0577
ode15i.m, ode15s.m: Transpose event outputs in BISTs (bug #59416).
Markus Mützel <markus.muetzel@gmx.de>
parents:
29041
diff
changeset
|
508 %! assert (sol.ie, [1, 2]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
509 %! assert (isfield (sol, "xe")); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
510 %! assert (isfield (sol, "ye")); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
511 %! 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
|
512 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
513 ## Events option, five output arguments |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
514 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
515 %! 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
|
516 %! [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
|
517 %! [-1e-4; 1e-4; 0], opt); |
29039
1aa69571a313
Use 1-based indexing for "ie" variable in ode15i/s (bug #59063).
Rik <rik@octave.org>
parents:
29025
diff
changeset
|
518 %! assert (t(end), 10, 1); |
1aa69571a313
Use 1-based indexing for "ie" variable in ode15i/s (bug #59063).
Rik <rik@octave.org>
parents:
29025
diff
changeset
|
519 %! assert (te, [10; 10], 0.2); |
1aa69571a313
Use 1-based indexing for "ie" variable in ode15i/s (bug #59063).
Rik <rik@octave.org>
parents:
29025
diff
changeset
|
520 %! assert (ie, [1; 2]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
521 |
24008
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
522 ## 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
|
523 %!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
|
524 %! 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
|
525 %! [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
|
526 %! [0, 1], [1, 1], [1, 1]); |
28929
9e43deb9bfc3
maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents:
28912
diff
changeset
|
527 %! assert (size (yout), [20, 2]); |
24008
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
528 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
529 %!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
|
530 %! 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
|
531 %! [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
|
532 %! [0, 1], [1, 1], [1; 1]); |
28929
9e43deb9bfc3
maint: Use semicolon after assert statement inside %!test blocks.
Rik <rik@octave.org>
parents:
28912
diff
changeset
|
533 %! assert (size (yout), [20, 2]); |
24008
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
534 |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
535 ## 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
|
536 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
537 %! 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
|
538 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
539 %! '"Jacobian" function must evaluate to a real square matrix'); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
540 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
541 ## 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
|
542 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
543 %! DFDY = [-0.04, 1; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
544 %! 0.04, 1]; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
545 %! DFDYP = [-1, 0, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
546 %! 0, -1, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
547 %! 0, 0, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
548 %! 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
|
549 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
550 %! '"Jacobian" matrices must be real square matrices'); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
551 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
552 ## 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
|
553 %!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
|
554 %! 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
|
555 %! 0.04, 1]); |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
556 %! 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
|
557 %! 0, -1, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
558 %! 0, 0, 0]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
559 %! 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
|
560 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
561 %! '"Jacobian" matrices must be real square matrices'); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
562 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
563 ## 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
|
564 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
565 %! 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
|
566 %! 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
|
567 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
568 %! '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
|
569 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
570 ## Jacobian single matrix |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
571 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
572 %! 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
|
573 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
574 %! '"Jacobian" field must be a function handle or 2-element cell array of square matrices'); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
575 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
576 ## 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
|
577 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
578 %! 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
|
579 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
580 %! '"Jacobian" field must be a function handle or 2-element cell array of square matrices'); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
581 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
582 ## Jacobian strange field |
28046
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
583 %!testif HAVE_SUNDIALS |
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
584 %! opt = odeset ("Jacobian", "_5yVNhWVJWJn47RKnzxPsyb_"); |
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
585 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
5664362da646
ode15i.m: Adapt BIST tests to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28045
diff
changeset
|
586 %! '"Jacobian" function "_5yVNhWVJWJn47RKnzxPsyb_" not found'); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
587 |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
588 %!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
|
589 %! ydot = [y - yp]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
590 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
591 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
592 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
593 %! 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
|
594 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
595 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
596 %! 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
|
597 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
598 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
599 %! 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
|
600 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
601 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
602 %! 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
|
603 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
604 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
605 %! 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
|
606 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
607 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
608 %! 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
|
609 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
610 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
611 %! 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
|
612 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
613 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
614 %! fail ("ode15i (@fun, 1, 1, 1)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
615 %! "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
|
616 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
617 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
618 %! 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
|
619 %! "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
|
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, 2], 1, [1, 2])", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
623 %! "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
|
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 |
25712
5625b2237a4d
* ode15i.m: In tests, use unlikely symbol name instead of "foo".
John W. Eaton <jwe@octave.org>
parents:
25113
diff
changeset
|
626 %! 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
|
627 %! 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
|
628 %! "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
|
629 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
630 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
631 %! 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
|
632 %! 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
|
633 %! "ode15i: RelTol must be scalar"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
634 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
635 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
636 %! 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
|
637 %! 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
|
638 %! "ode15i: RelTol must be positive"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
639 |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
640 %!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
|
641 %! opt = odeset ("AbsTol", "_5yVNhWVJWJn47RKnzxPsyb_"); |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
642 %! 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
|
643 %! "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
|
644 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
645 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
646 %! 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
|
647 %! 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
|
648 %! "ode15i: AbsTol must be positive"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
649 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
650 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
651 %! 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
|
652 %! fail ("[t, y] = ode15i (@fun, [0, 2], 2, 2, opt)", |
28047
bd7d82fc65ba
ode15i.m: Adapt another BIST test to new input validation in cset 13dba3c069f8.
Rik <rik@octave.org>
parents:
28046
diff
changeset
|
653 %! '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
|
654 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
655 %!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
|
656 %! 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
|
657 %! 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
|
658 %! "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
|
659 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
660 %!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
|
661 %! 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
|
662 %! 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
|
663 %! "ode15i: YP0 must be a numeric vector"); |