Mercurial > octave
annotate scripts/ode/ode15i.m @ 25712:5625b2237a4d stable
* ode15i.m: In tests, use unlikely symbol name instead of "foo".
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 01 Aug 2018 12:24:21 -0400 |
parents | 476fc012d09a |
children | 23483673ba43 |
rev | line source |
---|---|
25054
6652d3823428
maint: Update copyright dates in all source files.
John W. Eaton <jwe@octave.org>
parents:
25026
diff
changeset
|
1 ## Copyright (C) 2016-2018 Francesco Faccio <francesco.faccio@mail.polimi.it> |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
2 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
3 ## This file is part of Octave. |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
4 ## |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
24008
diff
changeset
|
5 ## 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
|
6 ## 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
|
7 ## 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
|
8 ## (at your option) any later version. |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
9 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
10 ## 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
|
11 ## 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
|
12 ## 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
|
13 ## GNU General Public License for more details. |
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 ## 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
|
16 ## 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
|
17 ## <https://www.gnu.org/licenses/>. |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
18 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
20 ## @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
|
21 ## @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
|
22 ## @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
|
23 ## @deftypefnx {} {@var{solution} =} ode15i (@dots{}) |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
24 ## @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
|
25 ## 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
|
26 ## index 1 Differential Algebraic Equations (DAEs). |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
27 ## |
25113
476fc012d09a
doc: Shorten very long first sentences of docstrings (bug #53388).
John W. Eaton <jwe@octave.org>
parents:
25054
diff
changeset
|
28 ## @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
|
29 ## 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
|
30 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
31 ## @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
|
32 ## 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
|
33 ## 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
|
34 ## 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
|
35 ## 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
|
36 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
37 ## @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
|
38 ## 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
|
39 ## 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
|
40 ## 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
|
41 ## instances. |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
42 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
43 ## @var{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
|
44 ## 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
|
45 ## 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
|
46 ## value in @var{y0} and @var{yp0}. |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
47 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
48 ## @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
|
49 ## @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
|
50 ## 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
|
51 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
52 ## 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
|
53 ## 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
|
54 ## |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
55 ## 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
|
56 ## 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
|
57 ## 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
|
58 ## 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
|
59 ## |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
60 ## 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
|
61 ## 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
|
62 ## 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
|
63 ## 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
|
64 ## @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
|
65 ## additional information returned. |
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 ## If no output arguments are requested, and no @code{OutputFcn} is specified |
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
68 ## in @var{ode_opt}, then the @code{OutputFcn} is set to @code{odeplot} and the |
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
69 ## 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
|
70 ## |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
71 ## 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
|
72 ## 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
|
73 ## @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
|
74 ## 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
|
75 ## of multiple Event functions. |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
76 ## |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
77 ## Example: Solve @nospell{Robertson's} equations: |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
78 ## |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
79 ## @smallexample |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
80 ## @group |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
81 ## 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
|
82 ## 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
|
83 ## -(@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
|
84 ## @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
|
85 ## endfunction |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
86 ## [@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
|
87 ## @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
|
88 ## @end smallexample |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
89 ## @seealso{decic, odeset, odeget} |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
90 ## @end deftypefn |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
91 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
92 function varargout = ode15i (fun, trange, y0, yp0, varargin) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
93 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
94 if (nargin < 4) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
95 print_usage (); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
96 endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
97 |
25026
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
98 solver = "ode15i"; |
f886561f9696
doc: improve differential eqtn docs and mention ode15i/ode15s (bug #51965).
Colin Macdonald <cbm@m.fsf.org>
parents:
24534
diff
changeset
|
99 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
100 n = numel (y0); |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
101 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
102 if (nargin > 4) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
103 options = varargin{1}; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
104 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
105 options = odeset (); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
106 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
107 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
108 ## Check fun, trange, y0, yp0 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
109 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
|
110 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
111 if (! isempty (options.Jacobian)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
112 if (ischar (options.Jacobian)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
113 try |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
114 options.Jacobian = str2func (options.Jacobian); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
115 catch |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
116 warning (lasterr); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
117 end_try_catch |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
118 if (! isa (options.Jacobian, "function_handle")) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
119 error ("Octave:invalid-input-arg", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
120 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
121 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
122 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
123 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
124 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
125 if (! isempty (options.OutputFcn)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
126 if (ischar (options.OutputFcn)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
127 try |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
128 options.OutputFcn = str2func (options.OutputFcn); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
129 catch |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
130 warning (lasterr); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
131 end_try_catch |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
132 if (! isa (options.OutputFcn, "function_handle")) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
133 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
134 [solver ": invalid value assigned to field 'OutputFcn'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
135 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
136 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
137 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
138 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
139 if (! isempty (options.Events)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
140 if (ischar (options.Events)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
141 try |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
142 options.Events = str2func (options.Events); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
143 catch |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
144 warning (lasterr); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
145 end_try_catch |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
146 if (! isa (options.Events, "function_handle") |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
147 && ! ismatrix (options.Events)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
148 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
149 [solver ": invalid value assigned to field 'Events'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
150 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
151 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
152 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
153 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
154 [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
|
155 |
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
|
156 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
|
157 "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
|
158 "MassSingular", "InitialSlope", "BDF"}; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
159 |
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
|
160 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
|
161 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
|
162 attributes = rmfield (attributes, ignorefields); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
163 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
164 classes = odeset (classes, "Vectorized", {}); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
165 attributes = odeset (attributes, "Jacobian", {}, "Vectorized", {}); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
166 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
167 options = odemergeopts ("ode15i", options, defaults, |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
168 classes, attributes, solver); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
169 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
170 ## Jacobian |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
171 options.havejac = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
172 options.havejacsparse = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
173 options.havejacfun = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
174 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
175 if (! isempty (options.Jacobian)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
176 options.havejac = true; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
177 if (iscell (options.Jacobian)) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
178 if (numel (options.Jacobian) == 2) |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
179 if (issparse (options.Jacobian{1}) && issparse (options.Jacobian{2})) |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
180 options.havejacsparse = true; # Jac is sparse cell |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
181 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
182 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
183 if (any (size (options.Jacobian{1}) != [n n]) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
184 || any (size (options.Jacobian{2}) != [n n]) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
185 || ! isnumeric (options.Jacobian{1}) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
186 || ! isnumeric (options.Jacobian{2}) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
187 || ! isreal (options.Jacobian{1}) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
188 || ! isreal (options.Jacobian{2})) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
189 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
190 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
191 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
192 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
193 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
194 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
195 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
196 |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
197 elseif (isa (options.Jacobian, "function_handle")) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
198 options.havejacfun = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
199 if (nargin (options.Jacobian) == 3) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
200 [A, B] = options.Jacobian (trange(1), y0, yp0); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
201 if (issparse (A) && issparse (B)) |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
202 options.havejacsparse = true; # Jac is sparse fun |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
203 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
204 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
205 if (any (size (A) != [n n]) || any (size (B) != [n n]) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
206 || ! isnumeric (A) || ! isnumeric (B) || ! isreal (A) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
207 || ! isreal (B)) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
208 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
209 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
210 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
211 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
212 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
213 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
214 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
215 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
216 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
217 [solver ": invalid value assigned to field 'Jacobian'"]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
218 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
219 endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
220 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
221 ## Abstol and Reltol |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
222 options.haveabstolvec = false; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
223 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
224 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
|
225 error ("Octave:invalid-input-arg", |
22941
d92ec2901770
Make ode15i,ode15s doc consistent with other ode functions.
Rik <rik@octave.org>
parents:
22938
diff
changeset
|
226 [solver ": invalid value assigned to field 'AbsTol'"]); |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
227 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
228 elseif (numel (options.AbsTol) == n) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
229 options.haveabstolvec = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
230 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
231 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
232 ## Stats |
22935
c9344df03da5
Allow case-insensitive option argument 'on' to ode solvers (bug #49918).
Rik <rik@octave.org>
parents:
22933
diff
changeset
|
233 options.havestats = strcmpi (options.Stats, "on"); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
234 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
235 ## 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
|
236 if (nargout == 1) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
237 options.Refine = 1; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
238 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
239 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
240 ## OutputFcn and OutputSel |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
241 if (isempty (options.OutputFcn) && nargout == 0) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
242 options.OutputFcn = @odeplot; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
243 options.haveoutputfunction = true; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
244 else |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
245 options.haveoutputfunction = ! isempty (options.OutputFcn); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
246 endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
247 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
248 options.haveoutputselection = ! isempty (options.OutputSel); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
249 if (options.haveoutputselection) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
250 options.OutputSel = options.OutputSel - 1; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
251 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
252 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
253 ## Events |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
254 options.haveeventfunction = ! isempty (options.Events); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
255 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
256 [t, y, te, ye, ie] = __ode15__ (fun, trange, y0, yp0, options); |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
257 |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
258 if (nargout == 2) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
259 varargout{1} = t; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
260 varargout{2} = y; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
261 elseif (nargout == 1) |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
262 varargout{1}.x = t; # Time stamps are saved in field x |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
263 varargout{1}.y = y; # Results are saved in field y |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
264 varargout{1}.solver = solver; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
265 if (options.haveeventfunction) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
266 varargout{1}.xe = te; # Time info when an event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
267 varargout{1}.ye = ye; # Results when an event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
268 varargout{1}.ie = ie; # Index info which event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
269 endif |
22938
a88ceac2aa53
Fix undefined return argument for more than 2 outputs from ode15i,ode15s.
Rik <rik@octave.org>
parents:
22935
diff
changeset
|
270 elseif (nargout > 2) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
271 varargout = cell (1,5); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
272 varargout{1} = t; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
273 varargout{2} = y; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
274 if (options.haveeventfunction) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
275 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
|
276 varargout{4} = ye; # Results when an event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
277 varargout{5} = ie; # Index info which event occurred |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
278 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
279 endif |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
280 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
281 endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
282 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
283 |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
284 %!demo |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
285 %! ## 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
|
286 %! 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
|
287 %! -(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
|
288 %! y(1) + y(2) + y(3) - 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
289 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
290 %! 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
|
291 %! y0 = [1; 0; 0]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
292 %! yp0 = [-1e-4; 1e-4; 0]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
293 %! 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
|
294 %! |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
295 %! [t, y] = ode15i (fun, tspan, y0, yp0, opt); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
296 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
297 %! y(:,2) = 1e4 * y(:, 2); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
298 %! figure (2); |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
299 %! semilogx (t, y, "o"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
300 %! xlabel ("time"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
301 %! ylabel ("species concentration"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
302 %! 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
|
303 %! legend ("y1", "y2", "y3"); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
304 |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
305 %!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
|
306 %! 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
|
307 %! -(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
|
308 %! y(1) + y(2) + y(3) - 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
309 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
310 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
311 %!function ref = fref () |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
312 %! 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
|
313 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
314 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
315 %!function ref2 = fref2 () |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
316 %! ref2 = [4e6 0 0 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
317 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
318 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
319 %!function [DFDY, DFDYP] = jacfundense (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
320 %! 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
|
321 %! 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
|
322 %! 1, 1, 1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
323 %! DFDYP = [-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
324 %! 0, -1, 0; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
325 %! 0, 0, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
326 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
327 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
328 %!function [DFDY, DFDYP] = jacfunsparse (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
329 %! 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
|
330 %! 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
|
331 %! 1, 1, 1]); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
332 %! DFDYP = sparse ([-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
333 %! 0, -1, 0; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
334 %! 0, 0, 0]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
335 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
336 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
337 %!function [DFDY, DFDYP] = jacwrong (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
338 %! DFDY = [-0.04, 1e4*y(3); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
339 %! 0.04, -1e4*y(3)-6e7*y(2)]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
340 %! DFDYP = [-1, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
341 %! 0, -1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
342 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
343 %! |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
344 %!function [DFDY, DFDYP, A] = jacwrong2 (t, y, yp) |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
345 %! DFDY = [-0.04, 1e4*y(3), 1e4*y(2); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
346 %! 0.04, -1e4*y(3)-6e7*y(2), -1e4*y(2); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
347 %! 1, 1, 1]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
348 %! DFDYP = [-1, 0, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
349 %! 0, -1, 0; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
350 %! 0, 0, 0]; |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
351 %! A = DFDY; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
352 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
353 %! |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
354 %!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
|
355 %! isterminal = [0, 1]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
356 %! if (t < 1e1) |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
357 %! val = [-1, -2]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
358 %! else |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
359 %! val = [1, 3]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
360 %! endif |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
361 %! |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
362 %! direction = [1, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
363 %!endfunction |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
364 |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
365 ## anonymous function instead of real function |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
366 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
367 %! ref = 0.049787079136413; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
368 %! ff = @(t, u, udot) udot + 3 * u; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
369 %! [t, y] = ode15i (ff, 0:1, 1, -3); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
370 %! 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
|
371 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
372 ## function passed as string |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
373 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
374 %! [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
|
375 %! 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
|
376 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
377 ## solve in intermidiate step |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
378 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
379 %! [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
|
380 %! 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
|
381 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
382 ## numel(trange) = 2 final value |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
383 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
384 %! [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
|
385 %! 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
|
386 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
387 ## With empty options |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
388 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
389 %! opt = odeset(); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
390 %! [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
|
391 %! [-1e-4; 1e-4; 0], opt); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
392 %! assert ([t(end), y(end,:)], fref2, 1e-3); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
393 %! opt = odeset(); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
394 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
395 ## Without options |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
396 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
397 %! [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
|
398 %! 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
|
399 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
400 ## InitialStep option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
401 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
402 %! 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
|
403 %! [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
|
404 %! 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
|
405 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
406 ## MaxStep option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
407 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
408 %! 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
|
409 %! [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
|
410 %! 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
|
411 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
412 ## AbsTol scalar option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
413 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
414 %! 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
|
415 %! [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
|
416 %! 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
|
417 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
418 ## AbsTol scalar and RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
419 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
420 %! 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
|
421 %! [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
|
422 %! 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
|
423 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
424 ## AbsTol vector option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
425 %!testif HAVE_SUNDIALS |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
426 %! 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
|
427 %! [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
|
428 %! 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
|
429 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
430 ## AbsTol vector and RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
431 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
432 %! 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
|
433 %! [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
|
434 %! 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
|
435 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
436 ## RelTol option |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
437 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
438 %! 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
|
439 %! [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
|
440 %! 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
|
441 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
442 ## Jacobian fun dense |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
443 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
444 %! 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
|
445 %! [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
|
446 %! 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
|
447 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
448 ## Jacobian fun dense as string |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
449 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
450 %! 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
|
451 %! [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
|
452 %! 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
|
453 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
454 ## Jacobian fun sparse |
23023
63a12df71848
avoid sparse jacobian tests for ode15i and ode15s if IDAKLU is missing
John W. Eaton <jwe@octave.org>
parents:
22941
diff
changeset
|
455 %!testif HAVE_SUNDIALS_IDAKLU |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
456 %! 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
|
457 %! [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
|
458 %! 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
|
459 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
460 ## 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
|
461 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
462 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
463 %! 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
|
464 %! [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
|
465 %! [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
|
466 %! 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
|
467 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
468 ## 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
|
469 #%!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
470 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
471 %! 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
|
472 %! 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
|
473 %! [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
|
474 %! [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
|
475 %! 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
|
476 %! 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
|
477 |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
478 ## 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
|
479 #%!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
480 %! YPref = [-0.001135972751027; -0.000000027483627; 0.001136000234654]; |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
481 %! 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
|
482 %! [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
|
483 %! [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
|
484 %! 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
|
485 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
486 ## Refine |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
487 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
488 %! 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
|
489 %! [t, y] = ode15i (@rob, [0, 100], [1; 0; 0], [-1e-4; 1e-4; 0]); |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
490 %! [t2, y2] = ode15i (@rob, [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
|
491 %! 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
|
492 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
493 ## Refine ignored if numel (trange) > 2 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
494 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
495 %! opt = odeset ("Refine", 3); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
496 %! [t, y] = ode15i (@rob, [0, 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
|
497 %! [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
|
498 %! assert (numel (t2), numel (t)); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
499 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
500 ## 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
|
501 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
502 %! 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
|
503 %! 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
|
504 %! assert (isfield (sol, "ie")); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
505 %! assert (sol.ie, [0;1]); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
506 %! assert (isfield (sol, "xe")); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
507 %! assert (isfield (sol, "ye")); |
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
508 %! 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
|
509 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
510 ## Events option, five output arguments |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
511 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
512 %! 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
|
513 %! [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
|
514 %! [-1e-4; 1e-4; 0], opt); |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
515 %! assert ([t(end), te', ie'], [10, 10, 10, 0, 1], [1, 0.2, 0.2, 0, 0]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
516 |
24008
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
517 ## 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
|
518 %!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
|
519 %! 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
|
520 %! [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
|
521 %! [0, 1], [1, 1], [1, 1]); |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
522 %! assert (size (yout), [20, 2]) |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
523 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
524 %!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
|
525 %! 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
|
526 %! [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
|
527 %! [0, 1], [1, 1], [1; 1]); |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
528 %! assert (size (yout), [20, 2]) |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
529 |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
530 ## 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
|
531 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
532 %! 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
|
533 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
534 %! "ode15i: invalid value assigned to field 'Jacobian'"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
535 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
536 ## 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
|
537 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
538 %! DFDY = [-0.04, 1; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
539 %! 0.04, 1]; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
540 %! DFDYP = [-1, 0, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
541 %! 0, -1, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
542 %! 0, 0, 0]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
543 %! 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
|
544 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
545 %! "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
|
546 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
547 ## Jacobian cell sparse wrong dimension |
23023
63a12df71848
avoid sparse jacobian tests for ode15i and ode15s if IDAKLU is missing
John W. Eaton <jwe@octave.org>
parents:
22941
diff
changeset
|
548 %!testif HAVE_SUNDIALS_IDAKLU |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
549 %! 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
|
550 %! 0.04, 1]); |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
551 %! 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
|
552 %! 0, -1, 0; |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
553 %! 0, 0, 0]); |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
554 %! 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
|
555 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
556 %! "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
|
557 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
558 ## 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
|
559 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
560 %! 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
|
561 %! 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
|
562 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
563 %! "invalid value assigned to field 'Jacobian'"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
564 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
565 ## Jacobian single matrix |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
566 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
567 %! 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
|
568 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
569 %! "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
|
570 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
571 ## 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
|
572 %!testif HAVE_SUNDIALS |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
573 %! 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
|
574 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
575 %! "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
|
576 |
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
577 ## Jacobian strange field |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
578 ## FIXME: we need a better way to silence the warning from odeset. |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
579 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
580 %! saved_opts = warning (); |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
581 %! warning ("off", "all"); |
25712
5625b2237a4d
* ode15i.m: In tests, use unlikely symbol name instead of "foo".
John W. Eaton <jwe@octave.org>
parents:
25113
diff
changeset
|
582 %! opt = odeset ("Jacobian", "_5yVNhWVJWJn47RKnzxPsyb_"); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
583 %! fail ("[t, y] = ode15i (@rob, [0, 4e6], [1; 0; 0], [-1e-4; 1e-4; 0], opt)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
584 %! "invalid value assigned to field 'Jacobian'"); |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
585 %! warning (saved_opts); |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
586 |
22910
23847979b91e
maint: Strip trailing whitespace from source files.
John W. Eaton <jwe@octave.org>
parents:
22902
diff
changeset
|
587 %!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
|
588 %! ydot = [y - yp]; |
22901
4c56f3ffec04
Add functions ode15i and ode15s
Francesco Faccio <francesco.faccio@mail.polimi.it>
parents:
diff
changeset
|
589 %!endfunction |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
590 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
591 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
592 %! 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
|
593 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
594 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
595 %! 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
|
596 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
597 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
598 %! 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
|
599 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
600 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
601 %! fail ("ode15i (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
|
602 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
603 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
604 %! fail ("ode15i (1, 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
|
605 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
606 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
607 %! fail ("ode15i (1, 1, 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
|
608 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
609 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
610 %! fail ("ode15i (1, 1, 1, 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
|
611 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
612 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
613 %! fail ("ode15i (@fun, 1, 1, 1)", |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
614 %! "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
|
615 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
616 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
617 %! 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
|
618 %! "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
|
619 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
620 %!testif HAVE_SUNDIALS |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
621 %! 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
|
622 %! "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
|
623 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
624 %!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
|
625 %! 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
|
626 %! 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
|
627 %! "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
|
628 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
629 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
630 %! 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
|
631 %! 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
|
632 %! "ode15i: RelTol must be scalar"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
633 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
634 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
635 %! 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
|
636 %! 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
|
637 %! "ode15i: RelTol must be positive"); |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
638 |
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
639 %!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
|
640 %! opt = odeset ("AbsTol", "_5yVNhWVJWJn47RKnzxPsyb_"); |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
641 %! 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
|
642 %! "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
|
643 |
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
644 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
645 %! 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
|
646 %! 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
|
647 %! "ode15i: AbsTol must be positive"); |
22928
dec22bceafa2
use testif HAVE_SUNDIALS for ode15 tests
John W. Eaton <jwe@octave.org>
parents:
22910
diff
changeset
|
648 |
22931
8133da976602
* ode15i.m: Make error tests conditional on HAVE_SUNDIALS.
John W. Eaton <jwe@octave.org>
parents:
22928
diff
changeset
|
649 %!testif HAVE_SUNDIALS |
22933
c3428bb9aca9
ode15i.m, ode15s.m: Follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22931
diff
changeset
|
650 %! 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
|
651 %! 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
|
652 %! "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
|
653 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
654 %!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
|
655 %! 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
|
656 %! 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
|
657 %! "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
|
658 |
6e7bb85e32b8
Allow ode15s and ode15i to work with row initial vectors (bug #50192).
Marco Caliari <marco.caliari@univr.it>
parents:
24007
diff
changeset
|
659 %!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
|
660 %! 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
|
661 %! 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
|
662 %! "ode15i: YP0 must be a numeric vector"); |