Mercurial > octave-nkf
annotate scripts/ode/ode45.m @ 20581:e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Use '!' instead of '~' for logical not.
Use '##' only for comments which standalone on a line.
Use double quotes instead of single quotes where possible.
Wrap long lines to < 80 characters.
* ode45.m, odeget.m, odeset.m: Use Octave coding conventions.
author | Rik <rik@octave.org> |
---|---|
date | Sat, 03 Oct 2015 22:10:45 -0700 |
parents | 25623ef2ff4f |
children | 45151de7423f |
rev | line source |
---|---|
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
1 ## Copyright (C) 2006-2012, Thomas Treichl <treichl@users.sourceforge.net> |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
2 ## Copyright (C) 2013, Roberto Porcu' <roberto.porcu@polimi.it> |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
3 ## Copyright (C) 2014, Jacopo Corno <jacopo.corno@gmail.com> |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
4 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
5 ## This file is part of Octave. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
6 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
7 ## Octave is free software; you can redistribute it and/or modify it |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
8 ## under the terms of the GNU General Public License as published by |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
9 ## the Free Software Foundation; either version 3 of the License, or (at |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
10 ## your option) any later version. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
11 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
12 ## Octave is distributed in the hope that it will be useful, but |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
13 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
14 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
15 ## General Public License for more details. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
16 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
17 ## You should have received a copy of the GNU General Public License |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
18 ## along with Octave; see the file COPYING. If not, see |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
19 ## <http://www.gnu.org/licenses/>. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
20 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
21 ## -*- texinfo -*- |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
22 ## @deftypefn {Function File} {[@var{t}, @var{y}] =} ode45 (@var{fun}, @var{trange}, @var{init}) |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
23 ## @deftypefnx {Function File} {[@var{t}, @var{y}] =} ode45 (@var{fun}, @var{trange}, @var{init}, @var{opt}) |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
24 ## @deftypefnx {Function File} {[@var{t}, @var{y}] =} ode45 (@dots{}, @var{par1}, @var{par2}, @dots{}) |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
25 ## @deftypefnx {Function File} {[@var{t}, @var{y}, @var{xe}, @var{ye}, @var{ie}] =} ode45 (@dots{}) |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
26 ## @deftypefnx {Function File} {@var{sol} =} ode45 (@var{fun}, @var{trange}, @var{init}, @dots{}) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
27 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
28 ## Solve a set of non-stiff Ordinary Differential Equations (non-stiff ODEs) |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
29 ## with the well known explicit Dormand-Prince method of order 4. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
30 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
31 ## The first input argument must be a function handle or inline function that |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
32 ## defines the ODE: @code{y' = f(t,y)}. The function must accept two inputs |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
33 ## where the first is time @var{t} and the second is a column vector of |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
34 ## unknowns @var{y}. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
35 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
36 ## @var{trange} specifies the time interval over which the ODE will be |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
37 ## evaluated. Usually, it is a two-element vector specifying the initial and |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
38 ## final times (@code{[tinit, tfinal]}). If there are more than two elements |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
39 ## then the solution will also be evaluated at these intermediate time |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
40 ## instances unless the integrate function called is |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
41 ## @command{integrate_n_steps}. If there is only one time value, then |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
42 ## @code{ode45} will raise an error unless the options structure has |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
43 ## non-empty fields named @var{"TimeStepNumber"} and @var{"TimeStepSize"}. |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
44 ## If the option @var{"TimeStepSize"} is not empty, then the stepper called |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
45 ## will be @command{integrate_const}. If @var{"TimeStepNumber"} is also |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
46 ## specified then the integrate function @command{integrate_n_steps} will be |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
47 ## used; otherwise, @command{integrate_adaptive} is used. For this last |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
48 ## possibility the user can set the tolerance for the timestep computation by |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
49 ## changing the option @var{"Tau"}, that has a default value of @math{1e-6}. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
50 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
51 ## The third input argument @var{init} contains the initial value for the |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
52 ## unknowns. If this is a row vector then the solution @var{y} will be a matrix |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
53 ## in which each column is the solution for the corresponding initial value |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
54 ## in @var{init}. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
55 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
56 ## If present, the fourth input argument specifies options to the ODE solver. |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
57 ## It is a structure typically generated by @code{odeset}. |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
58 ## |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
59 ## The function usually produces just two outputs. Variable @var{t} is a |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
60 ## column vector and contains the times where the solution was found. The |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
61 ## output @var{y} is a matrix in which each column refers to a different |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
62 ## unknown of the problem and each row corresponds to a time in @var{t}. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
63 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
64 ## For example, solve an anonymous implementation of the Van der Pol equation |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
65 ## |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
66 ## @example |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
67 ## @group |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
68 ## fvdp = @@(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)]; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
69 ## [T,Y] = ode45 (fvdp, [0 20], [2 0]); |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
70 ## @end group |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
71 ## @end example |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20575
diff
changeset
|
72 ## @seealso{odeset, odeget} |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
73 ## @end deftypefn |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
74 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
75 function varargout = ode45 (vfun, vslot, vinit, varargin) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
76 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
77 vorder = 5; # runge_kutta_45_dorpri uses local extrapolation |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
78 vsolver = "ode45"; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
79 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
80 if (nargin == 0) # Check number and types of all input arguments |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
81 help (vsolver); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
82 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
83 "number of input arguments must be greater than zero"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
84 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
85 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
86 if (nargin < 3) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
87 print_usage; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
88 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
89 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
90 if (nargin >= 4) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
91 if (! isstruct (varargin{1})) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
92 ## varargin{1:len} are parameters for vfun |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
93 vodeoptions = odeset; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
94 vodeoptions.vfunarguments = varargin; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
95 elseif (length (varargin) > 1) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
96 ## varargin{1} is an OdePkg options structure vopt |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
97 vodeoptions = odepkg_structure_check (varargin{1}, "ode45"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
98 vodeoptions.vfunarguments = {varargin{2:length(varargin)}}; |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
99 else # if (isstruct (varargin{1})) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
100 vodeoptions = odepkg_structure_check (varargin{1}, "ode45"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
101 vodeoptions.vfunarguments = {}; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
102 endif |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
103 else # nargin == 3 |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
104 vodeoptions = odeset; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
105 vodeoptions.vfunarguments = {}; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
106 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
107 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
108 if (! isvector (vslot) || ! isnumeric (vslot)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
109 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
110 "second input argument must be a valid vector"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
111 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
112 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
113 if (length (vslot) < 2 && ... |
20575
3339c9bdfe6a
Activate FSAL property in dorpri timestepper
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
114 (isempty (vodeoptions.TimeStepSize) ... |
3339c9bdfe6a
Activate FSAL property in dorpri timestepper
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
115 || isempty (vodeoptions.TimeStepNumber))) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
116 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
117 "second input argument must be a valid vector"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
118 elseif (vslot(2) == vslot(1)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
119 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
120 "second input argument must be a valid vector"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
121 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
122 vodeoptions.vdirection = sign (vslot(2) - vslot(1)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
123 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
124 vslot = vslot(:); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
125 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
126 if (! isvector (vinit) || ! isnumeric (vinit)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
127 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
128 "third input argument must be a valid numerical value"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
129 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
130 vinit = vinit(:); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
131 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
132 if (! (isa (vfun, "function_handle") || isa (vfun, "inline"))) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
133 error ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
134 "first input argument must be a valid function handle"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
135 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
136 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
137 ## Start preprocessing, have a look which options are set in vodeoptions, |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
138 ## check if an invalid or unused option is set |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
139 if (isempty (vodeoptions.TimeStepNumber) ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
140 && isempty (vodeoptions.TimeStepSize)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
141 integrate_func = "adaptive"; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
142 vodeoptions.vstepsizefixed = false; |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
143 elseif (! isempty (vodeoptions.TimeStepNumber) ... |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
144 && ! isempty (vodeoptions.TimeStepSize)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
145 integrate_func = "n_steps"; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
146 vodeoptions.vstepsizefixed = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
147 if (sign (vodeoptions.TimeStepSize) != vodeoptions.vdirection) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
148 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
149 "option 'TimeStepSize' has a wrong sign", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
150 "it will be corrected automatically"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
151 vodeoptions.TimeStepSize = (-1)*vodeoptions.TimeStepSize; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
152 endif |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
153 elseif (isempty (vodeoptions.TimeStepNumber) ... |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
154 && ! isempty (vodeoptions.TimeStepSize)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
155 integrate_func = "const"; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
156 vodeoptions.vstepsizefixed = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
157 if (sign (vodeoptions.TimeStepSize) != vodeoptions.vdirection) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
158 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
159 "option 'TimeStepSize' has a wrong sign", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
160 "it will be corrected automatically"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
161 vodeoptions.TimeStepSize = (-1)*vodeoptions.TimeStepSize; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
162 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
163 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
164 warning ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
165 "assuming an adaptive integrate function"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
166 integrate_func = "adaptive"; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
167 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
168 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
169 ## Get the default options that can be set with "odeset" temporarily |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
170 vodetemp = odeset; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
171 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
172 ## Implementation of the option RelTol has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
173 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
174 if (isempty (vodeoptions.RelTol) && ! vodeoptions.vstepsizefixed) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
175 vodeoptions.RelTol = 1e-3; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
176 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
177 "Option 'RelTol' not set, new value %f is used", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
178 vodeoptions.RelTol); |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
179 elseif (! isempty (vodeoptions.RelTol) && vodeoptions.vstepsizefixed) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
180 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
181 "Option 'RelTol' will be ignored if fixed time stamps are given"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
182 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
183 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
184 ## Implementation of the option AbsTol has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
185 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
186 if (isempty (vodeoptions.AbsTol) && ! vodeoptions.vstepsizefixed) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
187 vodeoptions.AbsTol = 1e-6; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
188 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
189 "Option 'AbsTol' not set, new value %f is used", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
190 vodeoptions.AbsTol); |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
191 elseif (! isempty (vodeoptions.AbsTol) && vodeoptions.vstepsizefixed) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
192 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
193 "Option 'AbsTol' will be ignored if fixed time stamps are given"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
194 else |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
195 vodeoptions.AbsTol = vodeoptions.AbsTol(:); # Create column vector |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
196 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
197 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
198 ## Implementation of the option NormControl has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
199 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
200 if (strcmp (vodeoptions.NormControl, "on")) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
201 vodeoptions.vnormcontrol = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
202 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
203 vodeoptions.vnormcontrol = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
204 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
205 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
206 ## Implementation of the option NonNegative has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
207 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
208 if (! isempty (vodeoptions.NonNegative)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
209 if (isempty (vodeoptions.Mass)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
210 vodeoptions.vhavenonnegative = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
211 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
212 vodeoptions.vhavenonnegative = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
213 warning ("OdePkg:InvalidArgument", ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
214 "Option 'NonNegative' will be ignored if mass matrix is set"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
215 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
216 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
217 vodeoptions.vhavenonnegative = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
218 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
219 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
220 ## Implementation of the option OutputFcn has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
221 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
222 if (isempty (vodeoptions.OutputFcn) && nargout == 0) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
223 vodeoptions.OutputFcn = @odeplot; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
224 vodeoptions.vhaveoutputfunction = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
225 elseif (isempty (vodeoptions.OutputFcn)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
226 vodeoptions.vhaveoutputfunction = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
227 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
228 vodeoptions.vhaveoutputfunction = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
229 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
230 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
231 ## Implementation of the option OutputSel has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
232 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
233 if (! isempty (vodeoptions.OutputSel)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
234 vodeoptions.vhaveoutputselection = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
235 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
236 vodeoptions.vhaveoutputselection = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
237 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
238 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
239 ## Implementation of the option OutputSave has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
240 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
241 if (isempty (vodeoptions.OutputSave)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
242 vodeoptions.OutputSave = 1; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
243 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
244 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
245 ## Implementation of the option Refine has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
246 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
247 if (vodeoptions.Refine > 0) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
248 vodeoptions.vhaverefine = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
249 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
250 vodeoptions.vhaverefine = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
251 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
252 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
253 ## Implementation of the option Stats has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
254 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
255 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
256 ## Implementation of the option InitialStep has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
257 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
258 if (isempty (vodeoptions.InitialStep) && strcmp (integrate_func, "adaptive")) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
259 vodeoptions.InitialStep = ... |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
260 vodeoptions.vdirection * starting_stepsize (vorder, vfun, vslot(1), vinit, |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
261 vodeoptions.AbsTol, |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
262 vodeoptions.RelTol, |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
263 vodeoptions.vnormcontrol); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
264 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
265 "option 'InitialStep' not set, estimated value %f is used", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
266 vodeoptions.InitialStep); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
267 elseif(isempty (vodeoptions.InitialStep)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
268 vodeoptions.InitialStep = odeget (vodeoptions, "TimeStepSize"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
269 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
270 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
271 ## Implementation of the option MaxStep has been finished. This option |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
272 ## can be set by the user to another value than default value. |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
273 if (isempty (vodeoptions.MaxStep) && ! vodeoptions.vstepsizefixed) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
274 vodeoptions.MaxStep = abs (vslot(end) - vslot(1)) / 10; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
275 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
276 "Option 'MaxStep' not set, new value %f is used", ... |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
277 vodeoptions.MaxStep); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
278 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
279 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
280 ## Implementation of the option Events has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
281 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
282 if (! isempty (vodeoptions.Events)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
283 vodeoptions.vhaveeventfunction = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
284 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
285 vodeoptions.vhaveeventfunction = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
286 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
287 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
288 ## The options "Jacobian", "JPattern" and "Vectorized" will be ignored |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
289 ## by this solver because this solver uses an explicit Runge-Kutta method |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
290 ## and therefore no Jacobian calculation is necessary. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
291 if (! isequal (vodeoptions.Jacobian, vodetemp.Jacobian)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
292 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
293 "option 'Jacobian' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
294 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
295 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
296 if (! isequal (vodeoptions.JPattern, vodetemp.JPattern)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
297 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
298 "option 'JPattern' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
299 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
300 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
301 if (! isequal (vodeoptions.Vectorized, vodetemp.Vectorized)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
302 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
303 "option 'Vectorized' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
304 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
305 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
306 if (! isequal (vodeoptions.NewtonTol, vodetemp.NewtonTol)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
307 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
308 "option 'NewtonTol' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
309 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
310 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
311 if (! isequal (vodeoptions.MaxNewtonIterations,vodetemp.MaxNewtonIterations)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
312 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
313 "option 'MaxNewtonIterations' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
314 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
315 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
316 ## Implementation of the option Mass has been finished. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
317 ## This option can be set by the user to another value than default value. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
318 if (! isempty (vodeoptions.Mass) && isnumeric (vodeoptions.Mass)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
319 vhavemasshandle = false; |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
320 vmass = vodeoptions.Mass; # constant mass |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
321 elseif (isa (vodeoptions.Mass, "function_handle")) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
322 vhavemasshandle = true; # mass defined by a function handle |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
323 else # no mass matrix - creating a diag-matrix of ones for mass |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
324 vhavemasshandle = false; # vmass = diag (ones (length (vinit), 1), 0); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
325 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
326 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
327 ## Implementation of the option MStateDependence has been finished. |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
328 ## This option can be set by the user to another value than default value. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
329 if (strcmp (vodeoptions.MStateDependence, "none")) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
330 vmassdependence = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
331 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
332 vmassdependence = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
333 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
334 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
335 ## Other options that are not used by this solver. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
336 ## Print a warning message to tell the user that the option is ignored. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
337 if (! isequal (vodeoptions.MvPattern, vodetemp.MvPattern)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
338 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
339 "option 'MvPattern' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
340 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
341 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
342 if (! isequal (vodeoptions.MassSingular, vodetemp.MassSingular)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
343 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
344 "option 'MassSingular' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
345 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
346 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
347 if (! isequal (vodeoptions.InitialSlope, vodetemp.InitialSlope)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
348 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
349 "option 'InitialSlope' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
350 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
351 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
352 if (! isequal (vodeoptions.MaxOrder, vodetemp.MaxOrder)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
353 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
354 "option 'MaxOrder' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
355 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
356 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
357 if (! isequal (vodeoptions.BDF, vodetemp.BDF)) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
358 warning ("OdePkg:InvalidArgument", ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
359 "option 'BDF' will be ignored by this solver"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
360 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
361 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
362 ## Starting the initialisation of the core solver ode45 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
363 SubOpts = vodeoptions; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
364 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
365 if (vhavemasshandle) # Handle only the dynamic mass matrix, |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
366 if (vmassdependence) # constant mass matrices have already |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
367 vmass = @(t,x) vodeoptions.Mass (t, x, vodeoptions.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
368 vfun = @(t,x) vmass (t, x, vodeoptions.vfunarguments{:}) ... |
20575
3339c9bdfe6a
Activate FSAL property in dorpri timestepper
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
369 \ vfun (t, x, vodeoptions.vfunarguments{:}); |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
370 else # if (vmassdependence == false) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
371 vmass = @(t) vodeoptions.Mass (t, vodeoptions.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
372 vfun = @(t,x) vmass (t, vodeoptions.vfunarguments{:}) ... |
20575
3339c9bdfe6a
Activate FSAL property in dorpri timestepper
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
373 \ vfun (t, x, vodeoptions.vfunarguments{:}); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
374 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
375 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
376 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
377 switch integrate_func |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
378 case "adaptive" |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
379 solution = integrate_adaptive (@runge_kutta_45_dorpri, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
380 vorder, vfun, vslot, vinit, SubOpts); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
381 case "n_steps" |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
382 solution = integrate_n_steps (@runge_kutta_45_dorpri, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
383 vfun, vslot(1), vinit, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
384 vodeoptions.TimeStepSize, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
385 vodeoptions.TimeStepNumber, SubOpts); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
386 case "const" |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
387 solution = integrate_const (@runge_kutta_45_dorpri, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
388 vfun, vslot, vinit, ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
389 vodeoptions.TimeStepSize, SubOpts); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
390 endswitch |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
391 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
392 ## Postprocessing, do whatever when terminating integration algorithm |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
393 if (vodeoptions.vhaveoutputfunction) # Cleanup plotter |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
394 feval (vodeoptions.OutputFcn, solution.t(end), ... |
20575
3339c9bdfe6a
Activate FSAL property in dorpri timestepper
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
395 solution.x(end,:)', "done", vodeoptions.vfunarguments{:}); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
396 endif |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
397 if (vodeoptions.vhaveeventfunction) # Cleanup event function handling |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
398 odepkg_event_handle (vodeoptions.Events, solution.t(end), ... |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
399 solution.x(end,:)', "done", |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
400 vodeoptions.vfunarguments{:}); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
401 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
402 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
403 ## Print additional information if option Stats is set |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
404 if (strcmp (vodeoptions.Stats, "on")) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
405 vhavestats = true; |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
406 vnsteps = solution.vcntloop-2; # vcntloop from 2..end |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
407 vnfailed = (solution.vcntcycles-1)-(solution.vcntloop-2)+1; # vcntcycl from 1..end |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
408 vnfevals = 7*(solution.vcntcycles-1); # number of ode evaluations |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
409 vndecomps = 0; # number of LU decompositions |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
410 vnpds = 0; # number of partial derivatives |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
411 vnlinsols = 0; # no. of linear systems solutions |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
412 ## Print cost statistics if no output argument is given |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
413 if (nargout == 0) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
414 vmsg = fprintf (1, "Number of successful steps: %d\n", vnsteps); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
415 vmsg = fprintf (1, "Number of failed attempts: %d\n", vnfailed); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
416 vmsg = fprintf (1, "Number of function calls: %d\n", vnfevals); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
417 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
418 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
419 vhavestats = false; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
420 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
421 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
422 if (nargout == 1) # Sort output variables, depends on nargout |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
423 varargout{1}.x = solution.t; # Time stamps are saved in field x |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
424 varargout{1}.y = solution.x; # Results are saved in field y |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
425 varargout{1}.solver = vsolver; # Solver name is saved in field solver |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
426 if (vodeoptions.vhaveeventfunction) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
427 varargout{1}.ie = solution.vevent{2}; # Index info which event occurred |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
428 varargout{1}.xe = solution.vevent{3}; # Time info when an event occurred |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
429 varargout{1}.ye = solution.vevent{4}; # Results when an event occurred |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
430 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
431 if (vhavestats) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
432 varargout{1}.stats = struct; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
433 varargout{1}.stats.nsteps = vnsteps; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
434 varargout{1}.stats.nfailed = vnfailed; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
435 varargout{1}.stats.nfevals = vnfevals; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
436 varargout{1}.stats.npds = vnpds; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
437 varargout{1}.stats.ndecomps = vndecomps; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
438 varargout{1}.stats.nlinsols = vnlinsols; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
439 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
440 elseif (nargout == 2) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
441 varargout{1} = solution.t; # Time stamps are first output argument |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
442 varargout{2} = solution.x; # Results are second output argument |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
443 elseif (nargout == 5) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
444 varargout{1} = solution.t; # Same as (nargout == 2) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
445 varargout{2} = solution.x; # Same as (nargout == 2) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
446 varargout{3} = []; # LabMat doesn't accept lines like |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
447 varargout{4} = []; # varargout{3} = varargout{4} = []; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
448 varargout{5} = []; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
449 if (vodeoptions.vhaveeventfunction) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
450 varargout{3} = solution.vevent{3}; # Time info when an event occurred |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
451 varargout{4} = solution.vevent{4}; # Results when an event occurred |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
452 varargout{5} = solution.vevent{2}; # Index info which event occurred |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
453 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
454 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
455 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
456 endfunction |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
457 |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
458 |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
459 ## We are using the "Van der Pol" implementation for all tests that are done |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
460 ## for this function. |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
461 ## For further tests we also define a reference solution (computed at high |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
462 ## accuracy) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
463 %!function [ydot] = fpol (vt, vy) # The Van der Pol |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
464 %! ydot = [vy(2); (1 - vy(1)^2) * vy(2) - vy(1)]; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
465 %!endfunction |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
466 %!function [vref] = fref () # The computed reference solution |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
467 %! vref = [0.32331666704577, -1.83297456798624]; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
468 %!endfunction |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
469 %!function [vjac] = fjac (vt, vy, varargin) # its Jacobian |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
470 %! vjac = [0, 1; -1 - 2 * vy(1) * vy(2), 1 - vy(1)^2]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
471 %!endfunction |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
472 %!function [vjac] = fjcc (vt, vy, varargin) # sparse type |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
473 %! vjac = sparse ([0, 1; -1 - 2 * vy(1) * vy(2), 1 - vy(1)^2]) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
474 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
475 %!function [vval, vtrm, vdir] = feve (vt, vy, varargin) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
476 %! vval = fpol (vt, vy, varargin); # We use the derivatives |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
477 %! vtrm = zeros (2,1); # that's why component 2 |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
478 %! vdir = ones (2,1); # seems to not be exact |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
479 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
480 %!function [vval, vtrm, vdir] = fevn (vt, vy, varargin) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
481 %! vval = fpol (vt, vy, varargin); # We use the derivatives |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
482 %! vtrm = ones (2,1); # that's why component 2 |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
483 %! vdir = ones (2,1); # seems to not be exact |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
484 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
485 %!function [vmas] = fmas (vt, vy, varargin) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
486 %! vmas = [1, 0; 0, 1]; # Dummy mass matrix for tests |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
487 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
488 %!function [vmas] = fmsa (vt, vy, varargin) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
489 %! vmas = sparse ([1, 0; 0, 1]); # A sparse dummy matrix |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
490 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
491 %!function [vout] = fout (vt, vy, vflag, varargin) |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
492 %! if (regexp (char (vflag), 'init') == 1) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
493 %! if (any (size (vt) != [2, 1])) error ('"fout" step "init"'); endif |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
494 %! elseif (isempty (vflag)) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
495 %! if (any (size (vt) != [1, 1])) error ('"fout" step "calc"'); endif |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
496 %! vout = false; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
497 %! elseif (regexp (char (vflag), 'done') == 1) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
498 %! if (any (size (vt) != [1, 1])) error ('"fout" step "done"'); endif |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
499 %! else |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
500 %! error ('"fout" invalid vflag'); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
501 %! endif |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
502 %!endfunction |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
503 %! |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
504 %! ## Turn off output of warning messages for all tests, turn them on |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
505 %! ## again if the last test is called |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
506 %!error # ouput argument |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
507 %! warning ("off", "OdePkg:InvalidArgument"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
508 %! B = ode45 (1, [0 25], [3 15 1]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
509 %!error # input argument number one |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
510 %! [vt, vy] = ode45 (1, [0 25], [3 15 1]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
511 %!error # input argument number two |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
512 %! [vt, vy] = ode45 (@fpol, 1, [3 15 1]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
513 %!test # two output arguments |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
514 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
515 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
516 %!test # not too many steps |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
517 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
518 %! assert (size (vt) < 20); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
519 %!test # anonymous function instead of real function |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
520 %! fvdb = @(vt,vy) [vy(2); (1 - vy(1)^2) * vy(2) - vy(1)]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
521 %! [vt, vy] = ode45 (fvdb, [0 2], [2 0]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
522 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
523 %!test # extra input arguments passed through |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
524 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], 12, 13, "KL"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
525 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
526 %!test # empty OdePkg structure *but* extra input arguments |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
527 %! vopt = odeset; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
528 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt, 12, 13, "KL"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
529 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
530 %!error # strange OdePkg structure |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
531 %! vopt = struct ("foo", 1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
532 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
533 %!test # Solve vdp in fixed step sizes |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
534 %! vopt = odeset("TimeStepSize", 0.1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
535 %! [vt, vy] = ode45 (@fpol, [0,2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
536 %! assert (vt(:), [0:0.1:2]', 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
537 %!test # Solve another anonymous function below zero |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
538 %! vref = [0, 14.77810590694212]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
539 %! [vt, vy] = ode45 (@(t,y) y, [-2 0], 2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
540 %! assert ([vt(end), vy(end,:)], vref, 1e-1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
541 %!test # InitialStep option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
542 %! vopt = odeset ("InitialStep", 1e-8); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
543 %! [vt, vy] = ode45 (@fpol, [0 0.2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
544 %! assert ([vt(2)-vt(1)], [1e-8], 1e-9); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
545 %!test # MaxStep option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
546 %! vopt = odeset ("MaxStep", 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
547 %! vsol = ode45 (@fpol, [0 0.2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
548 %! assert ([vsol.x(5)-vsol.x(4)], [1e-3], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
549 %!test # Solve with intermidiate step |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
550 %! vsol = ode45 (@fpol, [0 1 2], [2 0]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
551 %! assert (any((vsol.x-1) == 0)); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
552 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
553 %!test # Solve in backward direction starting at t=0 |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
554 %! vref = [-1.205364552835178, 0.951542399860817]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
555 %! vsol = ode45 (@fpol, [0 -2], [2 0]); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
556 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
557 %!test # Solve in backward direction starting at t=2 |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
558 %! vref = [-1.205364552835178, 0.951542399860817]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
559 %! vsol = ode45 (@fpol, [2 -2], fref); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
560 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
561 %!test # Solve in backward direction starting at t=2, with intermidiate step |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
562 %! vref = [-1.205364552835178, 0.951542399860817]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
563 %! vsol = ode45 (@fpol, [2 0 -2], fref); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
564 %! idx = find(vsol.x < 0, 1, "first") - 1; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
565 %! assert ([vsol.x(idx), vsol.y(idx,:)], [0 2 0], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
566 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
567 %!test # Solve another anonymous function in backward direction |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
568 %! vref = [-1, 0.367879437558975]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
569 %! vsol = ode45 (@(t,y) y, [0 -1], 1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
570 %! assert ([vsol.x(end), vsol.y(end,:)], vref, 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
571 %!test # Solve another anonymous function below zero |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
572 %! vref = [0, 14.77810590694212]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
573 %! vsol = ode45 (@(t,y) y, [-2 0], 2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
574 %! assert ([vsol.x(end), vsol.y(end,:)], vref, 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
575 %!test # Solve in backward direction starting at t=0 with MaxStep option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
576 %! vref = [-1.205364552835178, 0.951542399860817]; |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
577 %! vopt = odeset ("MaxStep", 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
578 %! vsol = ode45 (@fpol, [0 -2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
579 %! assert ([abs(vsol.x(8)-vsol.x(7))], [1e-3], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
580 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
581 %!test # AbsTol option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
582 %! vopt = odeset ("AbsTol", 1e-5); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
583 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
584 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
585 %!test # AbsTol and RelTol option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
586 %! vopt = odeset ("AbsTol", 1e-8, "RelTol", 1e-8); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
587 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
588 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
589 %!test # RelTol and NormControl option -- higher accuracy |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
590 %! vopt = odeset ("RelTol", 1e-8, "NormControl", "on"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
591 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
592 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-5); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
593 %!test # Keeps initial values while integrating |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
594 %! vopt = odeset ("NonNegative", 2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
595 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
596 %! assert ([vsol.x(end), vsol.y(end,:)], [2, 2, 0], 0.5); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
597 %!test # Details of OutputSel and Refine can't be tested |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
598 %! vopt = odeset ("OutputFcn", @fout, "OutputSel", 1, "Refine", 5); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
599 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
600 %!test # Details of OutputSave can't be tested |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
601 %! vopt = odeset ("OutputSave", 1, "OutputSel", 1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
602 %! vsla = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
603 %! vopt = odeset ("OutputSave", 2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
604 %! vslb = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
605 %! assert (length (vsla.x) + 1 >= 2 * length (vslb.x)) |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
606 %!test # Stats must add further elements in vsol |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
607 %! vopt = odeset ("Stats", "on"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
608 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
609 %! assert (isfield (vsol, "stats")); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
610 %! assert (isfield (vsol.stats, "nsteps")); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
611 %!test # Events option add further elements in vsol |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
612 %! vopt = odeset ("Events", @feve); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
613 %! vsol = ode45 (@fpol, [0 10], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
614 %! assert (isfield (vsol, "ie")); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
615 %! assert (vsol.ie(1), 2); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
616 %! assert (isfield (vsol, "xe")); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
617 %! assert (isfield (vsol, "ye")); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
618 %!test # Events option, now stop integration |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
619 %! vopt = odeset ("Events", @fevn, "NormControl", "on"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
620 %! vsol = ode45 (@fpol, [0 10], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
621 %! assert ([vsol.ie, vsol.xe, vsol.ye], ... |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
622 %! [2.0, 2.496110, -0.830550, -2.677589], 6e-1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
623 %!test # Events option, five output arguments |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
624 %! vopt = odeset ("Events", @fevn, "NormControl", "on"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
625 %! [vt, vy, vxe, vye, vie] = ode45 (@fpol, [0 10], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
626 %! assert ([vie, vxe, vye], ... |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
627 %! [2.0, 2.496110, -0.830550, -2.677589], 6e-1); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
628 %!test # Jacobian option |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
629 %! vopt = odeset ("Jacobian", @fjac); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
630 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
631 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
632 %!test # Jacobian option and sparse return value |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
633 %! vopt = odeset ("Jacobian", @fjcc); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
634 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
635 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
636 %! |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
637 %!## test for JPattern option is missing |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
638 %!## test for Vectorized option is missing |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
639 %!## test for NewtonTol option is missing |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
640 %!## test for MaxNewtonIterations option is missing |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
641 %! |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
642 %!test # Mass option as function |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
643 %! vopt = odeset ("Mass", @fmas); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
644 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
645 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
646 %!test # Mass option as matrix |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
647 %! vopt = odeset ("Mass", eye (2,2)); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
648 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
649 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
650 %!test # Mass option as sparse matrix |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
651 %! vopt = odeset ("Mass", sparse (eye (2,2))); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
652 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
653 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
654 %!test # Mass option as function and sparse matrix |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
655 %! vopt = odeset ("Mass", @fmsa); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
656 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
657 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
658 %!test # Mass option as function and MStateDependence |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
659 %! vopt = odeset ("Mass", @fmas, "MStateDependence", "strong"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
660 %! vsol = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
661 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
662 %!test # Set BDF option to something else than default |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
663 %! vopt = odeset ("BDF", "on"); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
664 %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt); |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
665 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-3); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
666 %! |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
667 %!## test for MvPattern option is missing |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
668 %!## test for InitialSlope option is missing |
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
669 %!## test for MaxOrder option is missing |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
670 %! |
20581
e368ce72a844
maint: Use Octave coding conventions for ode* functions.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
671 %! warning ("on", "OdePkg:InvalidArgument"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
672 |