]> Creatis software - CreaPhase.git/blob - octave_packages/odepkg-0.8.2/ode54d.m
Add a useful package (from Source forge) for octave
[CreaPhase.git] / octave_packages / odepkg-0.8.2 / ode54d.m
1 %# Copyright (C) 2008-2012, Thomas Treichl <treichl@users.sourceforge.net>
2 %# OdePkg - A package for solving ordinary differential equations and more
3 %#
4 %# This program is free software; you can redistribute it and/or modify
5 %# it under the terms of the GNU General Public License as published by
6 %# the Free Software Foundation; either version 2 of the License, or
7 %# (at your option) any later version.
8 %#
9 %# This program is distributed in the hope that it will be useful,
10 %# but WITHOUT ANY WARRANTY; without even the implied warranty of
11 %# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
12 %# GNU General Public License for more details.
13 %#
14 %# You should have received a copy of the GNU General Public License
15 %# along with this program; If not, see <http://www.gnu.org/licenses/>.
16
17 %# -*- texinfo -*-
18 %# @deftypefn  {Function File} {[@var{}] =} ode54d (@var{@@fun}, @var{slot}, @var{init}, @var{lags}, @var{hist}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
19 %# @deftypefnx {Command} {[@var{sol}] =} ode54d (@var{@@fun}, @var{slot}, @var{init}, @var{lags}, @var{hist}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
20 %# @deftypefnx {Command} {[@var{t}, @var{y}, [@var{xe}, @var{ye}, @var{ie}]] =} ode54d (@var{@@fun}, @var{slot}, @var{init}, @var{lags}, @var{hist}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
21 %#
22 %# This function file can be used to solve a set of non--stiff delay differential equations (non--stiff DDEs) with a modified version of the well known explicit Runge--Kutta method of order (2,3).
23 %#
24 %# If this function is called with no return argument then plot the solution over time in a figure window while solving the set of DDEs that are defined in a function and specified by the function handle @var{@@fun}. The second input argument @var{slot} is a double vector that defines the time slot, @var{init} is a double vector that defines the initial values of the states, @var{lags} is a double vector that describes the lags of time, @var{hist} is a double matrix and describes the history of the DDEs, @var{opt} can optionally be a structure array that keeps the options created with the command @command{odeset} and @var{par1}, @var{par2}, @dots{} can optionally be other input arguments of any type that have to be passed to the function defined by @var{@@fun}.
25 %#
26 %# In other words, this function will solve a problem of the form
27 %# @example
28 %# dy/dt = fun (t, y(t), y(t-lags(1), y(t-lags(2), @dots{})))
29 %# y(slot(1)) = init
30 %# y(slot(1)-lags(1)) = hist(1), y(slot(1)-lags(2)) = hist(2), @dots{} 
31 %# @end example
32 %#
33 %# If this function is called with one return argument then return the solution @var{sol} of type structure array after solving the set of DDEs. The solution @var{sol} has the fields @var{x} of type double column vector for the steps chosen by the solver, @var{y} of type double column vector for the solutions at each time step of @var{x}, @var{solver} of type string for the solver name and optionally the extended time stamp information @var{xe}, the extended solution information @var{ye} and the extended index information @var{ie} all of type double column vector that keep the informations of the event function if an event function handle is set in the option argument @var{opt}.
34 %#
35 %# If this function is called with more than one return argument then return the time stamps @var{t}, the solution values @var{y} and optionally the extended time stamp information @var{xe}, the extended solution information @var{ye} and the extended index information @var{ie} all of type double column vector.
36 %#
37 %# For example:
38 %# @itemize @minus
39 %# @item
40 %# the following code solves an anonymous implementation of a chaotic behavior
41 %#
42 %# @example
43 %# fcao = @@(vt, vy, vz) [2 * vz / (1 + vz^9.65) - vy];
44 %#
45 %# vopt = odeset ("NormControl", "on", "RelTol", 1e-3);
46 %# vsol = ode54d (fcao, [0, 100], 0.5, 2, 0.5, vopt);
47 %#
48 %# vlag = interp1 (vsol.x, vsol.y, vsol.x - 2);
49 %# plot (vsol.y, vlag); legend ("fcao (t,y,z)");
50 %# @end example
51 %#
52 %# @item
53 %# to solve the following problem with two delayed state variables
54 %#
55 %# @example
56 %# d y1(t)/dt = -y1(t)
57 %# d y2(t)/dt = -y2(t) + y1(t-5)
58 %# d y3(t)/dt = -y3(t) + y2(t-10)*y1(t-10)
59 %# @end example
60 %#
61 %# one might do the following
62 %#
63 %# @example
64 %# function f = fun (t, y, yd)
65 %# f(1) = -y(1);                   %% y1' = -y1(t)
66 %# f(2) = -y(2) + yd(1,1);         %% y2' = -y2(t) + y1(t-lags(1))
67 %# f(3) = -y(3) + yd(2,2)*yd(1,2); %% y3' = -y3(t) + y2(t-lags(2))*y1(t-lags(2))
68 %# endfunction
69 %# T = [0,20]
70 %# res = ode54d (@@fun, T, [1;1;1], [5, 10], ones (3,2));
71 %# @end example
72 %#
73 %# @end itemize
74 %# @end deftypefn
75 %#
76 %# @seealso{odepkg}
77
78 function [varargout] = ode54d (vfun, vslot, vinit, vlags, vhist, varargin)
79
80   if (nargin == 0) %# Check number and types of all input arguments
81     help ('ode54d');
82     error ('OdePkg:InvalidArgument', ...
83       'Number of input arguments must be greater than zero');
84
85   elseif (nargin < 5)
86     print_usage;
87
88   elseif (~isa (vfun, 'function_handle'))
89     error ('OdePkg:InvalidArgument', ...
90       'First input argument must be a valid function handle');
91
92   elseif (~isvector (vslot) || length (vslot) < 2)
93     error ('OdePkg:InvalidArgument', ...
94       'Second input argument must be a valid vector');
95
96   elseif (~isvector (vinit) || ~isnumeric (vinit))
97     error ('OdePkg:InvalidArgument', ...
98       'Third input argument must be a valid numerical value');
99
100   elseif (~isvector (vlags) || ~isnumeric (vlags))
101     error ('OdePkg:InvalidArgument', ...
102       'Fourth input argument must be a valid numerical value');
103
104   elseif ~(isnumeric (vhist) || isa (vhist, 'function_handle'))
105     error ('OdePkg:InvalidArgument', ...
106       'Fifth input argument must either be numeric or a function handle');
107
108   elseif (nargin >= 6)
109
110     if (~isstruct (varargin{1}))
111       %# varargin{1:len} are parameters for vfun
112       vodeoptions = odeset;
113       vfunarguments = varargin;
114
115     elseif (length (varargin) > 1)
116       %# varargin{1} is an OdePkg options structure vopt
117       vodeoptions = odepkg_structure_check (varargin{1}, 'ode54d');
118       vfunarguments = {varargin{2:length(varargin)}};
119
120     else %# if (isstruct (varargin{1}))
121       vodeoptions = odepkg_structure_check (varargin{1}, 'ode54d');
122       vfunarguments = {};
123
124     end
125
126   else %# if (nargin == 5)
127     vodeoptions = odeset; 
128     vfunarguments = {};
129   end
130
131   %# Start preprocessing, have a look which options have been set in
132   %# vodeoptions. Check if an invalid or unused option has been set and
133   %# print warnings.
134   vslot = vslot(:)'; %# Create a row vector
135   vinit = vinit(:)'; %# Create a row vector
136   vlags = vlags(:)'; %# Create a row vector
137
138   %# Check if the user has given fixed points of time
139   if (length (vslot) > 2), vstepsizegiven = true; %# Step size checking
140   else vstepsizegiven = false; end  
141
142   %# Get the default options that can be set with 'odeset' temporarily
143   vodetemp = odeset;
144
145   %# Implementation of the option RelTol has been finished. This option
146   %# can be set by the user to another value than default value.
147   if (isempty (vodeoptions.RelTol) && ~vstepsizegiven)
148     vodeoptions.RelTol = 1e-6;
149     warning ('OdePkg:InvalidOption', ...
150       'Option "RelTol" not set, new value %f is used', vodeoptions.RelTol);
151   elseif (~isempty (vodeoptions.RelTol) && vstepsizegiven)
152     warning ('OdePkg:InvalidOption', ...
153       'Option "RelTol" will be ignored if fixed time stamps are given');
154   %# This implementation has been added to odepkg_structure_check.m
155   %# elseif (~isscalar (vodeoptions.RelTol) && ~vstepsizegiven)
156   %# error ('OdePkg:InvalidOption', ...
157   %#   'Option "RelTol" must be set to a scalar value for this solver');
158   end
159
160   %# Implementation of the option AbsTol has been finished. This option
161   %# can be set by the user to another value than default value.
162   if (isempty (vodeoptions.AbsTol) && ~vstepsizegiven)
163     vodeoptions.AbsTol = 1e-6;
164     warning ('OdePkg:InvalidOption', ...
165       'Option "AbsTol" not set, new value %f is used', vodeoptions.AbsTol);
166   elseif (~isempty (vodeoptions.AbsTol) && vstepsizegiven)
167     warning ('OdePkg:InvalidOption', ...
168       'Option "AbsTol" will be ignored if fixed time stamps are given');
169   else %# create column vector
170     vodeoptions.AbsTol = vodeoptions.AbsTol(:);
171   end
172
173   %# Implementation of the option NormControl has been finished. This
174   %# option can be set by the user to another value than default value.
175   if (strcmp (vodeoptions.NormControl, 'on')), vnormcontrol = true;
176   else vnormcontrol = false;
177   end
178
179   %# Implementation of the option NonNegative has been finished. This
180   %# option can be set by the user to another value than default value.
181   if (~isempty (vodeoptions.NonNegative))
182     if (isempty (vodeoptions.Mass)), vhavenonnegative = true;
183     else
184       vhavenonnegative = false;
185       warning ('OdePkg:InvalidOption', ...
186         'Option "NonNegative" will be ignored if mass matrix is set');
187     end
188   else vhavenonnegative = false;
189   end
190
191   %# Implementation of the option OutputFcn has been finished. This
192   %# option can be set by the user to another value than default value.
193   if (isempty (vodeoptions.OutputFcn) && nargout == 0)
194     vodeoptions.OutputFcn = @odeplot;
195     vhaveoutputfunction = true;
196   elseif (isempty (vodeoptions.OutputFcn)), vhaveoutputfunction = false;
197   else vhaveoutputfunction = true;
198   end
199
200   %# Implementation of the option OutputSel has been finished. This
201   %# option can be set by the user to another value than default value.
202   if (~isempty (vodeoptions.OutputSel)), vhaveoutputselection = true;
203   else vhaveoutputselection = false; end
204
205   %# Implementation of the option Refine has been finished. This option
206   %# can be set by the user to another value than default value.
207   if (isequal (vodeoptions.Refine, vodetemp.Refine)), vhaverefine = true;
208   else vhaverefine = false; end
209
210   %# Implementation of the option Stats has been finished. This option
211   %# can be set by the user to another value than default value.
212
213   %# Implementation of the option InitialStep has been finished. This
214   %# option can be set by the user to another value than default value.
215   if (isempty (vodeoptions.InitialStep) && ~vstepsizegiven)
216     vodeoptions.InitialStep = abs (vslot(1,1) - vslot(1,2)) / 10;
217     vodeoptions.InitialStep = vodeoptions.InitialStep / 10^vodeoptions.Refine;
218     warning ('OdePkg:InvalidOption', ...
219       'Option "InitialStep" not set, new value %f is used', vodeoptions.InitialStep);
220   end
221
222   %# Implementation of the option MaxStep has been finished. This option
223   %# can be set by the user to another value than default value.
224   if (isempty (vodeoptions.MaxStep) && ~vstepsizegiven)
225     vodeoptions.MaxStep = abs (vslot(1,1) - vslot(1,length (vslot))) / 10;
226     %# vodeoptions.MaxStep = vodeoptions.MaxStep / 10^vodeoptions.Refine;
227     warning ('OdePkg:InvalidOption', ...
228       'Option "MaxStep" not set, new value %f is used', vodeoptions.MaxStep);
229   end
230
231   %# Implementation of the option Events has been finished. This option
232   %# can be set by the user to another value than default value.
233   if (~isempty (vodeoptions.Events)), vhaveeventfunction = true;
234   else vhaveeventfunction = false; end
235
236   %# The options 'Jacobian', 'JPattern' and 'Vectorized' will be ignored
237   %# by this solver because this solver uses an explicit Runge-Kutta
238   %# method and therefore no Jacobian calculation is necessary
239   if (~isequal (vodeoptions.Jacobian, vodetemp.Jacobian))
240     warning ('OdePkg:InvalidOption', ...
241       'Option "Jacobian" will be ignored by this solver');
242   end
243   if (~isequal (vodeoptions.JPattern, vodetemp.JPattern))
244     warning ('OdePkg:InvalidOption', ...
245       'Option "JPattern" will be ignored by this solver');
246   end
247   if (~isequal (vodeoptions.Vectorized, vodetemp.Vectorized))
248     warning ('OdePkg:InvalidOption', ...
249       'Option "Vectorized" will be ignored by this solver');
250   end
251   if (~isequal (vodeoptions.NewtonTol, vodetemp.NewtonTol))
252     warning ('OdePkg:InvalidArgument', ...
253       'Option "NewtonTol" will be ignored by this solver');
254   end
255   if (~isequal (vodeoptions.MaxNewtonIterations,...
256                 vodetemp.MaxNewtonIterations))
257     warning ('OdePkg:InvalidArgument', ...
258       'Option "MaxNewtonIterations" will be ignored by this solver');
259   end
260
261   %# Implementation of the option Mass has been finished. This option
262   %# can be set by the user to another value than default value.
263   if (~isempty (vodeoptions.Mass) && isnumeric (vodeoptions.Mass))
264     vhavemasshandle = false; vmass = vodeoptions.Mass; %# constant mass
265   elseif (isa (vodeoptions.Mass, 'function_handle'))
266     vhavemasshandle = true; %# mass defined by a function handle
267   else %# no mass matrix - creating a diag-matrix of ones for mass
268     vhavemasshandle = false; %# vmass = diag (ones (length (vinit), 1), 0);
269   end
270
271   %# Implementation of the option MStateDependence has been finished.
272   %# This option can be set by the user to another value than default
273   %# value. 
274   if (strcmp (vodeoptions.MStateDependence, 'none'))
275     vmassdependence = false;
276   else vmassdependence = true;
277   end
278
279   %# Other options that are not used by this solver. Print a warning
280   %# message to tell the user that the option(s) is/are ignored.
281   if (~isequal (vodeoptions.MvPattern, vodetemp.MvPattern))
282     warning ('OdePkg:InvalidOption', ...
283       'Option "MvPattern" will be ignored by this solver');
284   end
285   if (~isequal (vodeoptions.MassSingular, vodetemp.MassSingular))
286     warning ('OdePkg:InvalidOption', ...
287       'Option "MassSingular" will be ignored by this solver');
288   end
289   if (~isequal (vodeoptions.InitialSlope, vodetemp.InitialSlope))
290     warning ('OdePkg:InvalidOption', ...
291       'Option "InitialSlope" will be ignored by this solver');
292   end
293   if (~isequal (vodeoptions.MaxOrder, vodetemp.MaxOrder))
294     warning ('OdePkg:InvalidOption', ...
295       'Option "MaxOrder" will be ignored by this solver');
296   end
297   if (~isequal (vodeoptions.BDF, vodetemp.BDF))
298     warning ('OdePkg:InvalidOption', ...
299       'Option "BDF" will be ignored by this solver');
300   end
301
302   %# Starting the initialisation of the core solver ode54d 
303   vtimestamp  = vslot(1,1);           %# timestamp = start time
304   vtimelength = length (vslot);       %# length needed if fixed steps
305   vtimestop   = vslot(1,vtimelength); %# stop time = last value
306
307   if (~vstepsizegiven)
308     vstepsize = vodeoptions.InitialStep;
309     vminstepsize = (vtimestop - vtimestamp) / (1/eps);
310   else %# If step size is given then use the fixed time steps
311     vstepsize = abs (vslot(1,1) - vslot(1,2));
312     vminstepsize = eps; %# vslot(1,2) - vslot(1,1) - eps;
313   end
314
315   vretvaltime = vtimestamp; %# first timestamp output
316   if (vhaveoutputselection) %# first solution output
317     vretvalresult = vinit(vodeoptions.OutputSel);
318   else vretvalresult = vinit;
319   end
320
321   %# Initialize the OutputFcn
322   if (vhaveoutputfunction)
323     feval (vodeoptions.OutputFcn, vslot', ...
324       vretvalresult', 'init', vfunarguments{:});
325   end
326
327   %# Initialize the History
328   if (isnumeric (vhist))
329     vhmat = vhist;
330     vhavehistnumeric = true;
331   else %# it must be a function handle
332     for vcnt = 1:length (vlags);
333       vhmat(:,vcnt) = feval (vhist, (vslot(1)-vlags(vcnt)), vfunarguments{:});
334     end
335     vhavehistnumeric = false;
336   end
337
338   %# Initialize DDE variables for history calculation
339   vsaveddetime = [vtimestamp - vlags, vtimestamp]';
340   vsaveddeinput = [vhmat, vinit']';
341   vsavedderesult = [vhmat, vinit']';
342
343   %# Initialize the EventFcn
344   if (vhaveeventfunction)
345     odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
346       {vretvalresult', vhmat}, 'init', vfunarguments{:});
347   end
348
349   vpow = 1/5;                   %# 20071016, reported by Luis Randez
350   va = [0, 0, 0, 0, 0, 0;       %# The Dormand-Prince 5(4) coefficients
351         1/5, 0, 0, 0, 0, 0;     %# Coefficients proved on 20060827
352         3/40, 9/40, 0, 0, 0, 0; %# See p.91 in Ascher & Petzold
353         44/45, -56/15, 32/9, 0, 0, 0;
354         19372/6561, -25360/2187, 64448/6561, -212/729, 0, 0;
355         9017/3168, -355/33, 46732/5247, 49/176, -5103/18656, 0;
356         35/384, 0, 500/1113, 125/192, -2187/6784, 11/84];
357   %# 4th and 5th order b-coefficients
358   vb4 = [35/384; 0; 500/1113; 125/192; -2187/6784; 11/84; 0];
359   vb5 = [5179/57600; 0; 7571/16695; 393/640; -92097/339200; 187/2100; 1/40];
360   vc = sum (va, 2);
361
362   %# The solver main loop - stop if the endpoint has been reached
363   vcntloop = 2; vcntcycles = 1; vu = vinit; vk = vu' * zeros(1,7);
364   vcntiter = 0; vunhandledtermination = true;
365   while ((vtimestamp < vtimestop && vstepsize >= vminstepsize))
366
367     %# Hit the endpoint of the time slot exactely
368     if ((vtimestamp + vstepsize) > vtimestop)
369       vstepsize = vtimestop - vtimestamp; end
370
371     %# Estimate the seven results when using this solver
372     for j = 1:7
373       vthetime  = vtimestamp + vc(j,1) * vstepsize;
374       vtheinput = vu' + vstepsize * vk(:,1:j-1) * va(j,1:j-1)';
375       %# Claculate the history values (or get them from an external
376       %# function) that are needed for the next step of solving
377       if (vhavehistnumeric)
378         for vcnt = 1:length (vlags)
379           %# Direct implementation of a 'quadrature cubic Hermite interpolation'
380           %# found at the Faculty for Mathematics of the University of Stuttgart
381           %# http://mo.mathematik.uni-stuttgart.de/inhalt/aussage/aussage1269
382           vnumb = find (vthetime - vlags(vcnt) >= vsaveddetime);
383           velem = min (vnumb(end), length (vsaveddetime) - 1);
384           vstep = vsaveddetime(velem+1) - vsaveddetime(velem);
385           vdiff = (vthetime - vlags(vcnt) - vsaveddetime(velem)) / vstep;
386           vsubs = 1 - vdiff;
387           %# Calculation of the coefficients for the interpolation algorithm
388           vua = (1 + 2 * vdiff) * vsubs^2;
389           vub = (3 - 2 * vdiff) * vdiff^2;
390           vva = vstep * vdiff * vsubs^2;
391           vvb = -vstep * vsubs * vdiff^2;
392           vhmat(:,vcnt) = vua * vsaveddeinput(velem,:)' + ...
393               vub * vsaveddeinput(velem+1,:)' + ...
394               vva * vsavedderesult(velem,:)' + ...
395               vvb * vsavedderesult(velem+1,:)';
396         end
397       else %# the history must be a function handle
398         for vcnt = 1:length (vlags)
399           vhmat(:,vcnt) = feval ...
400             (vhist, vthetime - vlags(vcnt), vfunarguments{:});
401         end
402       end
403
404       if (vhavemasshandle)   %# Handle only the dynamic mass matrix,
405         if (vmassdependence) %# constant mass matrices have already
406           vmass = feval ...  %# been set before (if any)
407             (vodeoptions.Mass, vthetime, vtheinput, vfunarguments{:});
408         else                 %# if (vmassdependence == false)
409           vmass = feval ...  %# then we only have the time argument
410             (vodeoptions.Mass, vthetime, vfunarguments{:});
411         end
412         vk(:,j) = vmass \ feval ...
413           (vfun, vthetime, vtheinput, vhmat, vfunarguments{:});
414       else
415         vk(:,j) = feval ...
416           (vfun, vthetime, vtheinput, vhmat, vfunarguments{:});
417       end
418     end
419
420     %# Compute the 4th and the 5th order estimation
421     y4 = vu' + vstepsize * (vk * vb4);
422     y5 = vu' + vstepsize * (vk * vb5);
423     if (vhavenonnegative)
424       vu(vodeoptions.NonNegative) = abs (vu(vodeoptions.NonNegative));
425       y4(vodeoptions.NonNegative) = abs (y4(vodeoptions.NonNegative));
426       y5(vodeoptions.NonNegative) = abs (y5(vodeoptions.NonNegative));
427     end
428     vSaveVUForRefine = vu;
429
430     %# Calculate the absolute local truncation error and the acceptable error
431     if (~vstepsizegiven)
432       if (~vnormcontrol)
433         vdelta = y5 - y4;
434         vtau = max (vodeoptions.RelTol * vu', vodeoptions.AbsTol);
435       else
436         vdelta = norm (y5 - y4, Inf);
437         vtau = max (vodeoptions.RelTol * max (norm (vu', Inf), 1.0), ...
438                     vodeoptions.AbsTol);
439       end
440     else %# if (vstepsizegiven == true)
441       vdelta = 1; vtau = 2;
442     end
443
444     %# If the error is acceptable then update the vretval variables
445     if (all (vdelta <= vtau))
446       vtimestamp = vtimestamp + vstepsize;
447       vu = y5'; %# MC2001: the higher order estimation as "local extrapolation"
448       vretvaltime(vcntloop,:) = vtimestamp;
449       if (vhaveoutputselection)
450         vretvalresult(vcntloop,:) = vu(vodeoptions.OutputSel);
451       else
452         vretvalresult(vcntloop,:) = vu;
453       end
454       vcntloop = vcntloop + 1; vcntiter = 0;
455
456       %# Update DDE values for next history calculation      
457       vsaveddetime(end+1) = vtimestamp;
458       vsaveddeinput(end+1,:) = vtheinput';
459       vsavedderesult(end+1,:) = vu;
460
461       %# Call plot only if a valid result has been found, therefore this
462       %# code fragment has moved here. Stop integration if plot function
463       %# returns false
464       if (vhaveoutputfunction)
465         if (vhaverefine)                  %# Do interpolation
466           for vcnt = 0:vodeoptions.Refine %# Approximation between told and t
467             vapproxtime = (vcnt + 1) * vstepsize / (vodeoptions.Refine + 2);
468             vapproxvals = vSaveVUForRefine' + vapproxtime * (vk * vb5);
469             if (vhaveoutputselection)
470               vapproxvals = vapproxvals(vodeoptions.OutputSel);
471             end
472             feval (vodeoptions.OutputFcn, (vtimestamp - vstepsize) + vapproxtime, ...
473               vapproxvals, [], vfunarguments{:});
474           end
475         end
476         vpltret = feval (vodeoptions.OutputFcn, vtimestamp, ...
477           vretvalresult(vcntloop-1,:)', [], vfunarguments{:});
478         if (vpltret), vunhandledtermination = false; break; end
479       end
480
481       %# Call event only if a valid result has been found, therefore this
482       %# code fragment has moved here. Stop integration if veventbreak is
483       %# true
484       if (vhaveeventfunction)
485         vevent = ...
486           odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
487             {vu(:), vhmat}, [], vfunarguments{:});
488         if (~isempty (vevent{1}) && vevent{1} == 1)
489           vretvaltime(vcntloop-1,:) = vevent{3}(end,:);
490           vretvalresult(vcntloop-1,:) = vevent{4}(end,:);
491           vunhandledtermination = false; break;
492         end
493       end
494     end %# If the error is acceptable ...
495
496     %# Update the step size for the next integration step
497     if (~vstepsizegiven)
498       %# vdelta may be 0 or even negative - could be an iteration problem
499       vdelta = max (vdelta, eps); 
500       vstepsize = min (vodeoptions.MaxStep, ...
501         min (0.8 * vstepsize * (vtau ./ vdelta) .^ vpow));
502     elseif (vstepsizegiven)
503       if (vcntloop < vtimelength)
504         vstepsize = vslot(1,vcntloop-1) - vslot(1,vcntloop-2);
505       end
506     end
507
508     %# Update counters that count the number of iteration cycles
509     vcntcycles = vcntcycles + 1; %# Needed for postprocessing
510     vcntiter = vcntiter + 1;     %# Needed to find iteration problems
511
512     %# Stop solving because the last 1000 steps no successful valid
513     %# value has been found
514     if (vcntiter >= 5000)
515       error (['Solving has not been successful. The iterative', ...
516         ' integration loop exited at time t = %f before endpoint at', ...
517         ' tend = %f was reached. This happened because the iterative', ...
518         ' integration loop does not find a valid solution at this time', ...
519         ' stamp. Try to reduce the value of "InitialStep" and/or', ...
520         ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
521     end
522
523   end %# The main loop
524
525   %# Check if integration of the ode has been successful
526   if (vtimestamp < vtimestop)
527     if (vunhandledtermination == true)
528       error (['Solving has not been successful. The iterative', ...
529         ' integration loop exited at time t = %f', ...
530         ' before endpoint at tend = %f was reached. This may', ...
531         ' happen if the stepsize grows smaller than defined in', ...
532         ' vminstepsize. Try to reduce the value of "InitialStep" and/or', ...
533         ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
534     else
535       warning ('OdePkg:HideWarning', ...
536         ['Solver has been stopped by a call of "break" in', ...
537          ' the main iteration loop at time t = %f before endpoint at', ...
538          ' tend = %f was reached. This may happen because the @odeplot', ...
539          ' function returned "true" or the @event function returned "true".'], ...
540          vtimestamp, vtimestop);
541     end
542   end
543
544   %# Postprocessing, do whatever when terminating integration algorithm
545   if (vhaveoutputfunction) %# Cleanup plotter
546     feval (vodeoptions.OutputFcn, vtimestamp, ...
547       vretvalresult(vcntloop-1,:)', 'done', vfunarguments{:});
548   end
549   if (vhaveeventfunction)  %# Cleanup event function handling
550     odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
551       {vretvalresult(vcntloop-1,:), vhmat}, 'done', vfunarguments{:});
552   end
553
554   %# Print additional information if option Stats is set
555   if (strcmp (vodeoptions.Stats, 'on'))
556     vhavestats = true;
557     vnsteps    = vcntloop-2;                    %# vcntloop from 2..end
558     vnfailed   = (vcntcycles-1)-(vcntloop-2)+1; %# vcntcycl from 1..end
559     vnfevals   = 7*(vcntcycles-1);              %# number of ode evaluations
560     vndecomps  = 0;                             %# number of LU decompositions
561     vnpds      = 0;                             %# number of partial derivatives
562     vnlinsols  = 0;                             %# no. of solutions of linear systems
563     %# Print cost statistics if no output argument is given
564     if (nargout == 0)
565       vmsg = fprintf (1, 'Number of successful steps: %d', vnsteps);
566       vmsg = fprintf (1, 'Number of failed attempts:  %d', vnfailed);
567       vmsg = fprintf (1, 'Number of function calls:   %d', vnfevals);
568     end
569   else vhavestats = false;
570   end
571
572   if (nargout == 1)                 %# Sort output variables, depends on nargout
573     varargout{1}.x = vretvaltime;   %# Time stamps are saved in field x
574     varargout{1}.y = vretvalresult; %# Results are saved in field y
575     varargout{1}.solver = 'ode54d'; %# Solver name is saved in field solver
576     if (vhaveeventfunction) 
577       varargout{1}.ie = vevent{2};  %# Index info which event occured
578       varargout{1}.xe = vevent{3};  %# Time info when an event occured
579       varargout{1}.ye = vevent{4};  %# Results when an event occured
580     end
581     if (vhavestats)
582       varargout{1}.stats = struct;
583       varargout{1}.stats.nsteps   = vnsteps;
584       varargout{1}.stats.nfailed  = vnfailed;
585       varargout{1}.stats.nfevals  = vnfevals;
586       varargout{1}.stats.npds     = vnpds;
587       varargout{1}.stats.ndecomps = vndecomps;
588       varargout{1}.stats.nlinsols = vnlinsols;
589     end
590   elseif (nargout == 2)
591     varargout{1} = vretvaltime;     %# Time stamps are first output argument
592     varargout{2} = vretvalresult;   %# Results are second output argument
593   elseif (nargout == 5)
594     varargout{1} = vretvaltime;     %# Same as (nargout == 2)
595     varargout{2} = vretvalresult;   %# Same as (nargout == 2)
596     varargout{3} = [];              %# LabMat doesn't accept lines like
597     varargout{4} = [];              %# varargout{3} = varargout{4} = [];
598     varargout{5} = [];
599     if (vhaveeventfunction) 
600       varargout{3} = vevent{3};     %# Time info when an event occured
601       varargout{4} = vevent{4};     %# Results when an event occured
602       varargout{5} = vevent{2};     %# Index info which event occured
603     end
604   %# else nothing will be returned, varargout{1} undefined
605   end
606
607 %! # We are using a "pseudo-DDE" implementation for all tests that
608 %! # are done for this function. We also define an Events and a
609 %! # pseudo-Mass implementation. For further tests we also define a
610 %! # reference solution (computed at high accuracy) and an OutputFcn.
611 %!function [vyd] = fexp (vt, vy, vz, varargin)
612 %!  vyd(1,1) = exp (- vt) - vz(1); %# The DDEs that are
613 %!  vyd(2,1) = vy(1) - vz(2);      %# used for all examples
614 %!function [vval, vtrm, vdir] = feve (vt, vy, vz, varargin)
615 %!  vval = fexp (vt, vy, vz); %# We use the derivatives
616 %!  vtrm = zeros (2,1);       %# don't stop solving here
617 %!  vdir = ones (2,1);        %# in positive direction
618 %!function [vval, vtrm, vdir] = fevn (vt, vy, vz, varargin)
619 %!  vval = fexp (vt, vy, vz); %# We use the derivatives
620 %!  vtrm = ones (2,1);        %# stop solving here
621 %!  vdir = ones (2,1);        %# in positive direction
622 %!function [vmas] = fmas (vt, vy, vz, varargin)
623 %!  vmas =  [1, 0; 0, 1];     %# Dummy mass matrix for tests
624 %!function [vmas] = fmsa (vt, vy, vz, varargin)
625 %!  vmas = sparse ([1, 0; 0, 1]); %# A dummy sparse matrix
626 %!function [vref] = fref ()       %# The reference solution
627 %!  vref = [0.12194462133618, 0.01652432423938];
628 %!function [vout] = fout (vt, vy, vflag, varargin)
629 %!  if (regexp (char (vflag), 'init') == 1)
630 %!    if (any (size (vt) ~= [2, 1])) error ('"fout" step "init"'); end
631 %!  elseif (isempty (vflag))
632 %!    if (any (size (vt) ~= [1, 1])) error ('"fout" step "calc"'); end
633 %!    vout = false;
634 %!  elseif (regexp (char (vflag), 'done') == 1)
635 %!    if (any (size (vt) ~= [1, 1])) error ('"fout" step "done"'); end
636 %!  else error ('"fout" invalid vflag');
637 %!  end
638 %!
639 %! %# Turn off output of warning messages for all tests, turn them on
640 %! %# again if the last test is called
641 %!error %# input argument number one
642 %!  warning ('off', 'OdePkg:InvalidOption');
643 %!  B = ode54d (1, [0 5], [1; 0], 1, [1; 0]);
644 %!error %# input argument number two
645 %!  B = ode54d (@fexp, 1, [1; 0], 1, [1; 0]);
646 %!error %# input argument number three
647 %!  B = ode54d (@fexp, [0 5], 1, 1, [1; 0]);
648 %!error %# input argument number four
649 %!  B = ode54d (@fexp, [0 5], [1; 0], [1; 1], [1; 0]);
650 %!error %# input argument number five
651 %!  B = ode54d (@fexp, [0 5], [1; 0], 1, 1);
652 %!test %# one output argument
653 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0]);
654 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
655 %!  assert (isfield (vsol, 'solver'));
656 %!  assert (vsol.solver, 'ode54d');
657 %!test %# two output arguments
658 %!  [vt, vy] = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0]);
659 %!  assert ([vt(end), vy(end,:)], [5, fref], 1e-1);
660 %!test %# five output arguments and no Events
661 %!  [vt, vy, vxe, vye, vie] = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0]);
662 %!  assert ([vt(end), vy(end,:)], [5, fref], 1e-1);
663 %!  assert ([vie, vxe, vye], []);
664 %!test %# anonymous function instead of real function
665 %!  faym = @(vt, vy, vz) [exp(-vt) - vz(1); vy(1) - vz(2)];
666 %!  vsol = ode54d (faym, [0 5], [1; 0], 1, [1; 0]);
667 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
668 %!test %# extra input arguments passed trhough
669 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], 'KL');
670 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
671 %!test %# empty OdePkg structure *but* extra input arguments
672 %!  vopt = odeset;
673 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt, 12, 13, 'KL');
674 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
675 %!error %# strange OdePkg structure
676 %!  vopt = struct ('foo', 1);
677 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
678 %!test %# AbsTol option
679 %!  vopt = odeset ('AbsTol', 1e-5);
680 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
681 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
682 %!test %# AbsTol and RelTol option
683 %!  vopt = odeset ('AbsTol', 1e-7, 'RelTol', 1e-7);
684 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
685 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
686 %!test %# RelTol and NormControl option
687 %!  vopt = odeset ('AbsTol', 1e-7, 'NormControl', 'on');
688 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
689 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], .5e-1);
690 %!test %# NonNegative for second component
691 %!  vopt = odeset ('NonNegative', 1);
692 %!  vsol = ode54d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
693 %!  assert ([vsol.x(end), vsol.y(end,:)], [2.5, 0.001, 0.237], 1e-1);
694 %!test %# Details of OutputSel and Refine can't be tested
695 %!  vopt = odeset ('OutputFcn', @fout, 'OutputSel', 1, 'Refine', 5);
696 %!  vsol = ode54d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
697 %!test %# Stats must add further elements in vsol
698 %!  vopt = odeset ('Stats', 'on');
699 %!  vsol = ode54d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
700 %!  assert (isfield (vsol, 'stats'));
701 %!  assert (isfield (vsol.stats, 'nsteps'));
702 %!test %# InitialStep option
703 %!  vopt = odeset ('InitialStep', 1e-8);
704 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
705 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
706 %!test %# MaxStep option
707 %!  vopt = odeset ('MaxStep', 1e-2);
708 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
709 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
710 %!test %# Events option add further elements in vsol
711 %!  vopt = odeset ('Events', @feve);
712 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
713 %!  assert (isfield (vsol, 'ie'));
714 %!  assert (vsol.ie, [1; 1]);
715 %!  assert (isfield (vsol, 'xe'));
716 %!  assert (isfield (vsol, 'ye'));
717 %!test %# Events option, now stop integration
718 %!  warning ('off', 'OdePkg:HideWarning');
719 %!  vopt = odeset ('Events', @fevn, 'NormControl', 'on');
720 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
721 %!  assert ([vsol.ie, vsol.xe, vsol.ye], ...
722 %!    [1.0000, 2.9219, -0.2127, -0.2671], 1e-1);
723 %!test %# Events option, five output arguments
724 %!  vopt = odeset ('Events', @fevn, 'NormControl', 'on');
725 %!  [vt, vy, vxe, vye, vie] = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
726 %!  assert ([vie, vxe, vye], ...
727 %!    [1.0000, 2.9219, -0.2127, -0.2671], 1e-1);
728 %!
729 %! %# test for Jacobian option is missing
730 %! %# test for Jacobian (being a sparse matrix) is missing
731 %! %# test for JPattern option is missing
732 %! %# test for Vectorized option is missing
733 %! %# test for NewtonTol option is missing
734 %! %# test for MaxNewtonIterations option is missing
735 %!
736 %!test %# Mass option as function
737 %!  vopt = odeset ('Mass', eye (2,2));
738 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
739 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
740 %!test %# Mass option as matrix
741 %!  vopt = odeset ('Mass', eye (2,2));
742 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
743 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
744 %!test %# Mass option as sparse matrix
745 %!  vopt = odeset ('Mass', sparse (eye (2,2)));
746 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
747 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
748 %!test %# Mass option as function and sparse matrix
749 %!  vopt = odeset ('Mass', @fmsa);
750 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
751 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
752 %!test %# Mass option as function and MStateDependence
753 %!  vopt = odeset ('Mass', @fmas, 'MStateDependence', 'strong');
754 %!  vsol = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
755 %!  assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 1e-1);
756 %!test %# Set BDF option to something else than default
757 %!  vopt = odeset ('BDF', 'on');
758 %!  [vt, vy] = ode54d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
759 %!  assert ([vt(end), vy(end,:)], [5, fref], 0.5);
760 %!
761 %! %# test for MvPattern option is missing
762 %! %# test for InitialSlope option is missing
763 %! %# test for MaxOrder option is missing
764 %!
765 %!  warning ('on', 'OdePkg:InvalidOption');
766
767 %# Local Variables: ***
768 %# mode: octave ***
769 %# End: ***