annotate scripts/ode/ode15s.m @ 30564:796f54d4ddbf stable

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