Mercurial > octave-nkf
annotate scripts/ode/private/integrate_n_steps.m @ 20639:a260a6acb70f
fix test failures introduced by a22d8a2eb0e5
* scripts/ode/private/integrate_adaptive.m: fix stepping backwards, fix
invocation of OutputFcn, fix text of some error messages
* scripts/ode/private/integrate_const.m: remove use of option OutputSave
* scripts/ode/private/integrate_n_steps.m: remove use of option OutputSave
author | Carlo de Falco <carlo.defalco@polimi.it> |
---|---|
date | Sun, 11 Oct 2015 23:09:01 +0200 |
parents | b7ac1e94266e |
children |
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) 2013, Roberto Porcu' <roberto.porcu@polimi.it> |
20570
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
2 ## |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
3 ## This file is part of Octave. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
4 ## |
20570
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
5 ## Octave is free software; you can redistribute it and/or modify it |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
6 ## under the terms of the GNU General Public License as published by |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
7 ## the Free Software Foundation; either version 3 of the License, or (at |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
8 ## your option) any later version. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
9 ## |
20570
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
10 ## Octave is distributed in the hope that it will be useful, but |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
13 ## General Public License for more details. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
14 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
15 ## You should have received a copy of the GNU General Public License |
20570
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
16 ## along with Octave; see the file COPYING. If not, see |
6256f6e366ac
Fix copyright text in private ode functions
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20568
diff
changeset
|
17 ## <http://www.gnu.org/licenses/>. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
18 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
20 ## @deftypefn {Function File} {[@var{t}, @var{y}] =} integrate_n_steps (@var{@@stepper}, @var{@@fun}, @var{t0}, @var{x0}, @var{dt}, @var{n}, @var{options}) |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
21 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
22 ## This function file can be called by an ODE solver function in order to |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
23 ## integrate the set of ODEs on the interval @var{[t0,t0 + n*dt]} with a |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
24 ## constant timestep dt and on a fixed number of steps. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
25 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
26 ## This function must be called with two output arguments: @var{t} and @var{y}. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
27 ## Variable @var{t} is a column vector and contains the time stamps, instead |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
28 ## @var{y} is a matrix in which each column refers to a different unknown of |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
29 ## the problem and the rows number is the same of @var{t} rows number so that |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
30 ## each row of @var{y} contains the values of all unknowns at the time |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
31 ## value contained in the corresponding row in @var{t}. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
32 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
33 ## The first input argument must be a function_handle or an inline function |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
34 ## representing the stepper, that is the function responsible for step-by-step |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
35 ## integration. This function discriminates one method from the others. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
36 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
37 ## The second input argument is the order of the stepper. It is needed to |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
38 ## compute the adaptive timesteps. |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
39 ## |
20580
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
40 ## The third input argument is a function_handle or an inline function that |
25623ef2ff4f
doc: Rewrite docstrings for ode* family of functions.
Rik <rik@octave.org>
parents:
20570
diff
changeset
|
41 ## defines the set of ODE: |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
42 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
43 ## @ifhtml |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
44 ## @example |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
45 ## @math{y' = f(t,y)} |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
46 ## @end example |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
47 ## @end ifhtml |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
48 ## @ifnothtml |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
49 ## @math{y' = f(t,y)}. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
50 ## @end ifnothtml |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
51 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
52 ## The third input argument is the starting point for the integration. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
53 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
54 ## The fourth argument contains the initial conditions for the ODEs. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
55 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
56 ## The fifth input argument represents the fixed timestep while the sixth |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
57 ## contains the number of integration steps. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
58 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
59 ## The last argument is a struct with the options that may be |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
60 ## needed by the stepper. |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
61 ## @end deftypefn |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
62 ## |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
63 ## @seealso{integrate_adaptive, integrate_const} |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
64 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
65 function solution = integrate_n_steps (stepper, func, t0, x0, dt, n, options) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
66 |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
67 solution = struct (); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
68 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
69 ## first values for time and solution |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
70 x = x0(:); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
71 t = t0; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
72 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
73 vdirection = odeget (options, "vdirection", [], "fast"); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
74 if (sign (dt) != vdirection) |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
75 error ("OdePkg:InvalidArgument", "option 'InitialStep' has a wrong sign"); |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
76 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
77 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
78 comp = 0.0; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
79 tk = t0; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
80 options.comp = comp; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
81 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
82 ## Initialize the OutputFcn |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
83 if (options.vhaveoutputfunction) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
84 if (options.vhaveoutputselection) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
85 solution.vretout = x(options.OutputSel,end); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
86 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
87 solution.vretout = x; |
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 feval (options.OutputFcn, tspan, solution.vretout, "init", |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
90 options.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
91 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
92 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
93 ## Initialize the EventFcn |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
94 if (options.vhaveeventfunction) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
95 odepkg_event_handle (options.Events, t(end), x, "init", |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
96 options.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
97 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
98 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
99 solution.vcntloop = 2; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
100 solution.vcntcycles = 1; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
101 #vu = vinit; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
102 #vk = vu.' * zeros(1,6); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
103 vcntiter = 0; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
104 solution.vunhandledtermination = true; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
105 solution.vcntsave = 2; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
106 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
107 z = t; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
108 u = x; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
109 k_vals = feval (func, t , x, options.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
110 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
111 for i = 1:n |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
112 ## Compute the integration step from t to t+dt |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
113 [s, y, ~, k_vals] = stepper (func, z(end), u(:,end), dt, options, k_vals); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
114 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
115 [tk, comp] = kahan (tk, comp, dt); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
116 options.comp = comp; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
117 s(end) = tk; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
118 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
119 if (options.vhavenonnegative) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
120 x(options.NonNegative,end) = abs (x(options.NonNegative,end)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
121 y(options.NonNegative,end) = abs (y(options.NonNegative,end)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
122 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
123 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
124 if (options.vhaveoutputfunction && options.vhaverefine) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
125 vSaveVUForRefine = u(:,end); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
126 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
127 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
128 ## values on this interval for time and solution |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
129 z = [t(end);s]; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
130 u = [x(:,end),y]; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
131 |
20639
a260a6acb70f
fix test failures introduced by a22d8a2eb0e5
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20586
diff
changeset
|
132 x = [x,u(:,2:end)]; |
a260a6acb70f
fix test failures introduced by a22d8a2eb0e5
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20586
diff
changeset
|
133 t = [t;z(2:end)]; |
a260a6acb70f
fix test failures introduced by a22d8a2eb0e5
Carlo de Falco <carlo.defalco@polimi.it>
parents:
20586
diff
changeset
|
134 solution.vcntsave = solution.vcntsave + 1; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
135 solution.vcntloop = solution.vcntloop + 1; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
136 vcntiter = 0; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
137 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
138 ## Call plot only if a valid result has been found, therefore this code |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
139 ## fragment has moved here. Stop integration if plot function returns false |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
140 if (options.vhaveoutputfunction) |
20586
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
141 for vcnt = 0:options.Refine # Approximation between told and t |
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
142 if (options.vhaverefine) # Do interpolation |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
143 vapproxtime = (vcnt + 1) / (options.Refine + 2); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
144 vapproxvals = (1 - vapproxtime) * vSaveVUForRefine ... |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
145 + (vapproxtime) * y(:,end); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
146 vapproxtime = s(end) + vapproxtime*dt; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
147 else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
148 vapproxvals = x(:,end); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
149 vapproxtime = t(end); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
150 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
151 if (options.vhaveoutputselection) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
152 vapproxvals = vapproxvals(options.OutputSel); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
153 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
154 vpltret = feval (options.OutputFcn, vapproxtime, vapproxvals, [], |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
155 options.vfunarguments{:}); |
20586
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
156 if (vpltret) # Leave refinement loop |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
157 break; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
158 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
159 endfor |
20586
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
160 if (vpltret) # Leave main loop |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
161 solution.vunhandledtermination = false; |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
162 break; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
163 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
164 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
165 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
166 ## Call event only if a valid result has been found, therefore this |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
167 ## code fragment has moved here. Stop integration if veventbreak is |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
168 ## true |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
169 if (options.vhaveeventfunction) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
170 solution.vevent = odepkg_event_handle (options.Events, t(end), x(:,end), |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
171 [], options.vfunarguments{:}); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
172 if (! isempty (solution.vevent{1}) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
173 && solution.vevent{1} == 1) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
174 t(solution.vcntloop-1,:) = solution.vevent{3}(end,:); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
175 x(:,solution.vcntloop-1) = solution.vevent{4}(end,:)'; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
176 solution.vunhandledtermination = false; |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
177 break; |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
178 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
179 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
180 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
181 ## Update counters that count the number of iteration cycles |
20586
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
182 solution.vcntcycles = solution.vcntcycles + 1; # Needed for cost statistics |
b7ac1e94266e
maint: Further clean up of functions in ode/private dir.
Rik <rik@octave.org>
parents:
20584
diff
changeset
|
183 vcntiter = vcntiter + 1; # Needed to find iteration problems |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
184 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
185 ## Stop solving because the last 1000 steps no successful valid |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
186 ## value has been found |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
187 if (vcntiter >= 5000) |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
188 error (["Solving has not been successful. The iterative", |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
189 " integration loop exited at time t = %f before endpoint at", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
190 " tend = %f was reached. This happened because the iterative", |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
191 " integration loop does not find a valid solution at this time", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
192 " stamp. Try to reduce the value of 'InitialStep' and/or", |
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
193 " 'MaxStep' with the command 'odeset'.\n"], |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
194 s(end), tspan(end)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
195 endif |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
196 endfor |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
197 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
198 ## Check if integration of the ode has been successful |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
199 #if (vdirection * z(end) < vdirection * tspan(end)) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
200 # if (solution.vunhandledtermination == true) |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
201 # error ("OdePkg:InvalidArgument", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
202 # ["Solving has not been successful. The iterative", |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
203 # " integration loop exited at time t = %f", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
204 # " before endpoint at tend = %f was reached. This may", |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
205 # " happen if the stepsize grows smaller than defined in", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
206 # " vminstepsize. Try to reduce the value of 'InitialStep'", |
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
207 # " and/or 'MaxStep' with the command 'odeset'.\n"], |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
208 # z(end), tspan(end)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
209 # else |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
210 # warning ("OdePkg:InvalidArgument", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
211 # ["Solver has been stopped by a call of 'break' in the main", |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
212 # " iteration loop at time t = %f before endpoint at tend = %f", |
20584
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
213 # " was reached. This may happen because the @odeplot function", |
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
214 # " returned 'true' or the @event function returned", |
eb9e2d187ed2
maint: Use Octave coding conventions in scripts/ode/private dir.
Rik <rik@octave.org>
parents:
20580
diff
changeset
|
215 # " 'true'.\n"], |
20568
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
216 # z(end), tspan(end)); |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
217 # endif |
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 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
220 solution.t = t; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
221 solution.x = x'; |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
222 |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
223 endfunction |
fcb792acab9b
Moving ode45, odeset, odeget, and levenshtein from odepkg to core.
jcorno <jacopo.corno@gmail.com>
parents:
diff
changeset
|
224 |