1 %# Copyright (C) 2008-2012, Thomas Treichl <treichl@users.sourceforge.net>
2 %# OdePkg - A package for solving ordinary differential equations and more
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.
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.
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/>.
18 %# @deftypefn {Function File} {[@var{}] =} ode45d (@var{@@fun}, @var{slot}, @var{init}, @var{lags}, @var{hist}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
19 %# @deftypefnx {Command} {[@var{sol}] =} ode45d (@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}]] =} ode45d (@var{@@fun}, @var{slot}, @var{init}, @var{lags}, @var{hist}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
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 (4,5).
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}.
26 %# In other words, this function will solve a problem of the form
28 %# dy/dt = fun (t, y(t), y(t-lags(1), y(t-lags(2), @dots{})))
30 %# y(slot(1)-lags(1)) = hist(1), y(slot(1)-lags(2)) = hist(2), @dots{}
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}.
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.
40 %# the following code solves an anonymous implementation of a chaotic behavior
43 %# fcao = @@(vt, vy, vz) [2 * vz / (1 + vz^9.65) - vy];
45 %# vopt = odeset ("NormControl", "on", "RelTol", 1e-3);
46 %# vsol = ode45d (fcao, [0, 100], 0.5, 2, 0.5, vopt);
48 %# vlag = interp1 (vsol.x, vsol.y, vsol.x - 2);
49 %# plot (vsol.y, vlag); legend ("fcao (t,y,z)");
53 %# to solve the following problem with two delayed state variables
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)
61 %# one might do the following
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))
70 %# res = ode45d (@@fun, T, [1;1;1], [5, 10], ones (3,2));
78 function [varargout] = ode45d (vfun, vslot, vinit, vlags, vhist, varargin)
80 if (nargin == 0) %# Check number and types of all input arguments
82 error ('OdePkg:InvalidArgument', ...
83 'Number of input arguments must be greater than zero');
88 elseif (~isa (vfun, 'function_handle'))
89 error ('OdePkg:InvalidArgument', ...
90 'First input argument must be a valid function handle');
92 elseif (~isvector (vslot) || length (vslot) < 2)
93 error ('OdePkg:InvalidArgument', ...
94 'Second input argument must be a valid vector');
96 elseif (~isvector (vinit) || ~isnumeric (vinit))
97 error ('OdePkg:InvalidArgument', ...
98 'Third input argument must be a valid numerical value');
100 elseif (~isvector (vlags) || ~isnumeric (vlags))
101 error ('OdePkg:InvalidArgument', ...
102 'Fourth input argument must be a valid numerical value');
104 elseif ~(isnumeric (vhist) || isa (vhist, 'function_handle'))
105 error ('OdePkg:InvalidArgument', ...
106 'Fifth input argument must either be numeric or a function handle');
110 if (~isstruct (varargin{1}))
111 %# varargin{1:len} are parameters for vfun
112 vodeoptions = odeset;
113 vfunarguments = varargin;
115 elseif (length (varargin) > 1)
116 %# varargin{1} is an OdePkg options structure vopt
117 vodeoptions = odepkg_structure_check (varargin{1}, 'ode45d');
118 vfunarguments = {varargin{2:length(varargin)}};
120 else %# if (isstruct (varargin{1}))
121 vodeoptions = odepkg_structure_check (varargin{1}, 'ode45d');
126 else %# if (nargin == 5)
127 vodeoptions = odeset;
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
134 vslot = vslot(:)'; %# Create a row vector
135 vinit = vinit(:)'; %# Create a row vector
136 vlags = vlags(:)'; %# Create a row vector
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
142 %# Get the default options that can be set with 'odeset' temporarily
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');
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(:);
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;
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;
184 vhavenonnegative = false;
185 warning ('OdePkg:InvalidOption', ...
186 'Option "NonNegative" will be ignored if mass matrix is set');
188 else vhavenonnegative = false;
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;
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
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
210 %# Implementation of the option Stats has been finished. This option
211 %# can be set by the user to another value than default value.
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);
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);
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
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');
243 if (~isequal (vodeoptions.JPattern, vodetemp.JPattern))
244 warning ('OdePkg:InvalidOption', ...
245 'Option "JPattern" will be ignored by this solver');
247 if (~isequal (vodeoptions.Vectorized, vodetemp.Vectorized))
248 warning ('OdePkg:InvalidOption', ...
249 'Option "Vectorized" will be ignored by this solver');
251 if (~isequal (vodeoptions.NewtonTol, vodetemp.NewtonTol))
252 warning ('OdePkg:InvalidArgument', ...
253 'Option "NewtonTol" will be ignored by this solver');
255 if (~isequal (vodeoptions.MaxNewtonIterations,...
256 vodetemp.MaxNewtonIterations))
257 warning ('OdePkg:InvalidArgument', ...
258 'Option "MaxNewtonIterations" will be ignored by this solver');
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);
271 %# Implementation of the option MStateDependence has been finished.
272 %# This option can be set by the user to another value than default
274 if (strcmp (vodeoptions.MStateDependence, 'none'))
275 vmassdependence = false;
276 else vmassdependence = true;
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');
285 if (~isequal (vodeoptions.MassSingular, vodetemp.MassSingular))
286 warning ('OdePkg:InvalidOption', ...
287 'Option "MassSingular" will be ignored by this solver');
289 if (~isequal (vodeoptions.InitialSlope, vodetemp.InitialSlope))
290 warning ('OdePkg:InvalidOption', ...
291 'Option "InitialSlope" will be ignored by this solver');
293 if (~isequal (vodeoptions.MaxOrder, vodetemp.MaxOrder))
294 warning ('OdePkg:InvalidOption', ...
295 'Option "MaxOrder" will be ignored by this solver');
297 if (~isequal (vodeoptions.BDF, vodetemp.BDF))
298 warning ('OdePkg:InvalidOption', ...
299 'Option "BDF" will be ignored by this solver');
302 %# Starting the initialisation of the core solver ode45d
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
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;
315 vretvaltime = vtimestamp; %# first timestamp output
316 if (vhaveoutputselection) %# first solution output
317 vretvalresult = vinit(vodeoptions.OutputSel);
318 else vretvalresult = vinit;
321 %# Initialize the OutputFcn
322 if (vhaveoutputfunction)
323 feval (vodeoptions.OutputFcn, vslot', ...
324 vretvalresult', 'init', vfunarguments{:});
327 %# Initialize the History
328 if (isnumeric (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{:});
335 vhavehistnumeric = false;
338 %# Initialize DDE variables for history calculation
339 vsaveddetime = [vtimestamp - vlags, vtimestamp]';
340 vsaveddeinput = [vhmat, vinit']';
341 vsavedderesult = [vhmat, vinit']';
343 %# Initialize the EventFcn
344 if (vhaveeventfunction)
345 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
346 {vretvalresult', vhmat}, 'init', vfunarguments{:});
349 vpow = 1/5; %# 20071016, reported by Luis Randez
350 va = [0, 0, 0, 0, 0; %# The Runge-Kutta-Fehlberg 4(5) coefficients
351 1/4, 0, 0, 0, 0; %# Coefficients proved on 20060827
352 3/32, 9/32, 0, 0, 0; %# See p.91 in Ascher & Petzold
353 1932/2197, -7200/2197, 7296/2197, 0, 0;
354 439/216, -8, 3680/513, -845/4104, 0;
355 -8/27, 2, -3544/2565, 1859/4104, -11/40];
356 %# 4th and 5th order b-coefficients
357 vb4 = [25/216; 0; 1408/2565; 2197/4104; -1/5; 0];
358 vb5 = [16/135; 0; 6656/12825; 28561/56430; -9/50; 2/55];
361 %# The solver main loop - stop if endpoint has been reached
362 vcntloop = 2; vcntcycles = 1; vu = vinit; vk = vu' * zeros(1,6);
363 vcntiter = 0; vunhandledtermination = true;
364 while ((vtimestamp < vtimestop && vstepsize >= vminstepsize))
366 %# Hit the endpoint of the time slot exactely
367 if ((vtimestamp + vstepsize) > vtimestop)
368 vstepsize = vtimestop - vtimestamp; end
370 %# Estimate the six results when using this solver
372 vthetime = vtimestamp + vc(j,1) * vstepsize;
373 vtheinput = vu' + vstepsize * vk(:,1:j-1) * va(j,1:j-1)';
374 %# Claculate the history values (or get them from an external
375 %# function) that are needed for the next step of solving
376 if (vhavehistnumeric)
377 for vcnt = 1:length (vlags)
378 %# Direct implementation of a 'quadrature cubic Hermite interpolation'
379 %# found at the Faculty for Mathematics of the University of Stuttgart
380 %# http://mo.mathematik.uni-stuttgart.de/inhalt/aussage/aussage1269
381 vnumb = find (vthetime - vlags(vcnt) >= vsaveddetime);
382 velem = min (vnumb(end), length (vsaveddetime) - 1);
383 vstep = vsaveddetime(velem+1) - vsaveddetime(velem);
384 vdiff = (vthetime - vlags(vcnt) - vsaveddetime(velem)) / vstep;
386 %# Calculation of the coefficients for the interpolation algorithm
387 vua = (1 + 2 * vdiff) * vsubs^2;
388 vub = (3 - 2 * vdiff) * vdiff^2;
389 vva = vstep * vdiff * vsubs^2;
390 vvb = -vstep * vsubs * vdiff^2;
391 vhmat(:,vcnt) = vua * vsaveddeinput(velem,:)' + ...
392 vub * vsaveddeinput(velem+1,:)' + ...
393 vva * vsavedderesult(velem,:)' + ...
394 vvb * vsavedderesult(velem+1,:)';
396 else %# the history must be a function handle
397 for vcnt = 1:length (vlags)
398 vhmat(:,vcnt) = feval ...
399 (vhist, vthetime - vlags(vcnt), vfunarguments{:});
403 if (vhavemasshandle) %# Handle only the dynamic mass matrix,
404 if (vmassdependence) %# constant mass matrices have already
405 vmass = feval ... %# been set before (if any)
406 (vodeoptions.Mass, vthetime, vtheinput, vfunarguments{:});
407 else %# if (vmassdependence == false)
408 vmass = feval ... %# then we only have the time argument
409 (vodeoptions.Mass, vthetime, vfunarguments{:});
411 vk(:,j) = vmass \ feval ...
412 (vfun, vthetime, vtheinput, vhmat, vfunarguments{:});
415 (vfun, vthetime, vtheinput, vhmat, vfunarguments{:});
419 %# Compute the 4th and the 5th order estimation
420 y4 = vu' + vstepsize * (vk * vb4);
421 y5 = vu' + vstepsize * (vk * vb5);
422 if (vhavenonnegative)
423 vu(vodeoptions.NonNegative) = abs (vu(vodeoptions.NonNegative));
424 y4(vodeoptions.NonNegative) = abs (y4(vodeoptions.NonNegative));
425 y5(vodeoptions.NonNegative) = abs (y5(vodeoptions.NonNegative));
427 vSaveVUForRefine = vu;
429 %# Calculate the absolute local truncation error and the acceptable error
433 vtau = max (vodeoptions.RelTol * vu', vodeoptions.AbsTol);
435 vdelta = norm (y5 - y4, Inf);
436 vtau = max (vodeoptions.RelTol * max (norm (vu', Inf), 1.0), ...
439 else %# if (vstepsizegiven == true)
440 vdelta = 1; vtau = 2;
443 %# If the error is acceptable then update the vretval variables
444 if (all (vdelta <= vtau))
445 vtimestamp = vtimestamp + vstepsize;
446 vu = y5'; %# MC2001: the higher order estimation as "local extrapolation"
447 vretvaltime(vcntloop,:) = vtimestamp;
448 if (vhaveoutputselection)
449 vretvalresult(vcntloop,:) = vu(vodeoptions.OutputSel);
451 vretvalresult(vcntloop,:) = vu;
453 vcntloop = vcntloop + 1; vcntiter = 0;
455 %# Update DDE values for next history calculation
456 vsaveddetime(end+1) = vtimestamp;
457 vsaveddeinput(end+1,:) = vtheinput';
458 vsavedderesult(end+1,:) = vu;
460 %# Call plot only if a valid result has been found, therefore this
461 %# code fragment has moved here. Stop integration if plot function
463 if (vhaveoutputfunction)
464 if (vhaverefine) %# Do interpolation
465 for vcnt = 0:vodeoptions.Refine %# Approximation between told and t
466 vapproxtime = (vcnt + 1) * vstepsize / (vodeoptions.Refine + 2);
467 vapproxvals = vSaveVUForRefine' + vapproxtime * (vk * vb5);
468 if (vhaveoutputselection)
469 vapproxvals = vapproxvals(vodeoptions.OutputSel);
471 feval (vodeoptions.OutputFcn, (vtimestamp - vstepsize) + vapproxtime, ...
472 vapproxvals, [], vfunarguments{:});
475 vpltret = feval (vodeoptions.OutputFcn, vtimestamp, ...
476 vretvalresult(vcntloop-1,:)', [], vfunarguments{:});
477 if (vpltret), vunhandledtermination = false; break; end
480 %# Call event only if a valid result has been found, therefore this
481 %# code fragment has moved here. Stop integration if veventbreak is
483 if (vhaveeventfunction)
485 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
486 {vu(:), vhmat}, [], vfunarguments{:});
487 if (~isempty (vevent{1}) && vevent{1} == 1)
488 vretvaltime(vcntloop-1,:) = vevent{3}(end,:);
489 vretvalresult(vcntloop-1,:) = vevent{4}(end,:);
490 vunhandledtermination = false; break;
493 end %# If the error is acceptable ...
495 %# Update the step size for the next integration step
497 %# vdelta may be 0 or even negative - could be an iteration problem
498 vdelta = max (vdelta, eps);
499 vstepsize = min (vodeoptions.MaxStep, ...
500 min (0.8 * vstepsize * (vtau ./ vdelta) .^ vpow));
501 elseif (vstepsizegiven)
502 if (vcntloop < vtimelength)
503 vstepsize = vslot(1,vcntloop-1) - vslot(1,vcntloop-2);
507 %# Update counters that count the number of iteration cycles
508 vcntcycles = vcntcycles + 1; %# Needed for postprocessing
509 vcntiter = vcntiter + 1; %# Needed to find iteration problems
511 %# Stop solving because the last 1000 steps no successful valid
512 %# value has been found
513 if (vcntiter >= 5000)
514 error (['Solving has not been successful. The iterative', ...
515 ' integration loop exited at time t = %f before endpoint at', ...
516 ' tend = %f was reached. This happened because the iterative', ...
517 ' integration loop does not find a valid solution at this time', ...
518 ' stamp. Try to reduce the value of "InitialStep" and/or', ...
519 ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
524 %# Check if integration of the ode has been successful
525 if (vtimestamp < vtimestop)
526 if (vunhandledtermination == true)
527 error (['Solving has not been successful. The iterative', ...
528 ' integration loop exited at time t = %f', ...
529 ' before endpoint at tend = %f was reached. This may', ...
530 ' happen if the stepsize grows smaller than defined in', ...
531 ' vminstepsize. Try to reduce the value of "InitialStep" and/or', ...
532 ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
534 warning ('OdePkg:HideWarning', ...
535 ['Solver has been stopped by a call of "break" in', ...
536 ' the main iteration loop at time t = %f before endpoint at', ...
537 ' tend = %f was reached. This may happen because the @odeplot', ...
538 ' function returned "true" or the @event function returned "true".'], ...
539 vtimestamp, vtimestop);
543 %# Postprocessing, do whatever when terminating integration algorithm
544 if (vhaveoutputfunction) %# Cleanup plotter
545 feval (vodeoptions.OutputFcn, vtimestamp, ...
546 vretvalresult(vcntloop-1,:)', 'done', vfunarguments{:});
548 if (vhaveeventfunction) %# Cleanup event function handling
549 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
550 {vretvalresult(vcntloop-1,:), vhmat}, 'done', vfunarguments{:});
553 %# Print additional information if option Stats is set
554 if (strcmp (vodeoptions.Stats, 'on'))
556 vnsteps = vcntloop-2; %# vcntloop from 2..end
557 vnfailed = (vcntcycles-1)-(vcntloop-2)+1; %# vcntcycl from 1..end
558 vnfevals = 6*(vcntcycles-1); %# number of ode evaluations
559 vndecomps = 0; %# number of LU decompositions
560 vnpds = 0; %# number of partial derivatives
561 vnlinsols = 0; %# no. of solutions of linear systems
562 %# Print cost statistics if no output argument is given
564 vmsg = fprintf (1, 'Number of successful steps: %d', vnsteps);
565 vmsg = fprintf (1, 'Number of failed attempts: %d', vnfailed);
566 vmsg = fprintf (1, 'Number of function calls: %d', vnfevals);
568 else vhavestats = false;
571 if (nargout == 1) %# Sort output variables, depends on nargout
572 varargout{1}.x = vretvaltime; %# Time stamps are saved in field x
573 varargout{1}.y = vretvalresult; %# Results are saved in field y
574 varargout{1}.solver = 'ode45d'; %# Solver name is saved in field solver
575 if (vhaveeventfunction)
576 varargout{1}.ie = vevent{2}; %# Index info which event occured
577 varargout{1}.xe = vevent{3}; %# Time info when an event occured
578 varargout{1}.ye = vevent{4}; %# Results when an event occured
581 varargout{1}.stats = struct;
582 varargout{1}.stats.nsteps = vnsteps;
583 varargout{1}.stats.nfailed = vnfailed;
584 varargout{1}.stats.nfevals = vnfevals;
585 varargout{1}.stats.npds = vnpds;
586 varargout{1}.stats.ndecomps = vndecomps;
587 varargout{1}.stats.nlinsols = vnlinsols;
589 elseif (nargout == 2)
590 varargout{1} = vretvaltime; %# Time stamps are first output argument
591 varargout{2} = vretvalresult; %# Results are second output argument
592 elseif (nargout == 5)
593 varargout{1} = vretvaltime; %# Same as (nargout == 2)
594 varargout{2} = vretvalresult; %# Same as (nargout == 2)
595 varargout{3} = []; %# LabMat doesn't accept lines like
596 varargout{4} = []; %# varargout{3} = varargout{4} = [];
598 if (vhaveeventfunction)
599 varargout{3} = vevent{3}; %# Time info when an event occured
600 varargout{4} = vevent{4}; %# Results when an event occured
601 varargout{5} = vevent{2}; %# Index info which event occured
603 %# else nothing will be returned, varargout{1} undefined
606 %! # We are using a "pseudo-DDE" implementation for all tests that
607 %! # are done for this function. We also define an Events and a
608 %! # pseudo-Mass implementation. For further tests we also define a
609 %! # reference solution (computed at high accuracy) and an OutputFcn.
610 %!function [vyd] = fexp (vt, vy, vz, varargin)
611 %! vyd(1,1) = exp (- vt) - vz(1); %# The DDEs that are
612 %! vyd(2,1) = vy(1) - vz(2); %# used for all examples
613 %!function [vval, vtrm, vdir] = feve (vt, vy, vz, varargin)
614 %! vval = fexp (vt, vy, vz); %# We use the derivatives
615 %! vtrm = zeros (2,1); %# don't stop solving here
616 %! vdir = ones (2,1); %# in positive direction
617 %!function [vval, vtrm, vdir] = fevn (vt, vy, vz, varargin)
618 %! vval = fexp (vt, vy, vz); %# We use the derivatives
619 %! vtrm = ones (2,1); %# stop solving here
620 %! vdir = ones (2,1); %# in positive direction
621 %!function [vmas] = fmas (vt, vy, vz, varargin)
622 %! vmas = [1, 0; 0, 1]; %# Dummy mass matrix for tests
623 %!function [vmas] = fmsa (vt, vy, vz, varargin)
624 %! vmas = sparse ([1, 0; 0, 1]); %# A dummy sparse matrix
625 %!function [vref] = fref () %# The reference solution
626 %! vref = [0.12194462133618, 0.01652432423938];
627 %!function [vout] = fout (vt, vy, vflag, varargin)
628 %! if (regexp (char (vflag), 'init') == 1)
629 %! if (any (size (vt) ~= [2, 1])) error ('"fout" step "init"'); end
630 %! elseif (isempty (vflag))
631 %! if (any (size (vt) ~= [1, 1])) error ('"fout" step "calc"'); end
633 %! elseif (regexp (char (vflag), 'done') == 1)
634 %! if (any (size (vt) ~= [1, 1])) error ('"fout" step "done"'); end
635 %! else error ('"fout" invalid vflag');
638 %! %# Turn off output of warning messages for all tests, turn them on
639 %! %# again if the last test is called
640 %!error %# input argument number one
641 %! warning ('off', 'OdePkg:InvalidOption');
642 %! B = ode45d (1, [0 5], [1; 0], 1, [1; 0]);
643 %!error %# input argument number two
644 %! B = ode45d (@fexp, 1, [1; 0], 1, [1; 0]);
645 %!error %# input argument number three
646 %! B = ode45d (@fexp, [0 5], 1, 1, [1; 0]);
647 %!error %# input argument number four
648 %! B = ode45d (@fexp, [0 5], [1; 0], [1; 1], [1; 0]);
649 %!error %# input argument number five
650 %! B = ode45d (@fexp, [0 5], [1; 0], 1, 1);
651 %!test %# one output argument
652 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0]);
653 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
654 %! assert (isfield (vsol, 'solver'));
655 %! assert (vsol.solver, 'ode45d');
656 %!test %# two output arguments
657 %! [vt, vy] = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0]);
658 %! assert ([vt(end), vy(end,:)], [5, fref], 0.5);
659 %!test %# five output arguments and no Events
660 %! [vt, vy, vxe, vye, vie] = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0]);
661 %! assert ([vt(end), vy(end,:)], [5, fref], 0.5);
662 %! assert ([vie, vxe, vye], []);
663 %!test %# anonymous function instead of real function
664 %! faym = @(vt, vy, vz) [exp(-vt) - vz(1); vy(1) - vz(2)];
665 %! vsol = ode45d (faym, [0 5], [1; 0], 1, [1; 0]);
666 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
667 %!test %# extra input arguments passed trhough
668 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], 'KL');
669 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
670 %!test %# empty OdePkg structure *but* extra input arguments
672 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt, 12, 13, 'KL');
673 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
674 %!error %# strange OdePkg structure
675 %! vopt = struct ('foo', 1);
676 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
677 %!test %# AbsTol option
678 %! vopt = odeset ('AbsTol', 1e-5);
679 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
680 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
681 %!test %# AbsTol and RelTol option
682 %! vopt = odeset ('AbsTol', 1e-7, 'RelTol', 1e-7);
683 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
684 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
685 %!test %# RelTol and NormControl option
686 %! vopt = odeset ('AbsTol', 1e-7, 'NormControl', 'on');
687 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
688 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], .5e-1);
689 %!test %# NonNegative for second component
690 %! vopt = odeset ('NonNegative', 1);
691 %! vsol = ode45d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
692 %! assert ([vsol.x(end), vsol.y(end,:)], [2.5, 0.001, 0.237], 0.5);
693 %!test %# Details of OutputSel and Refine can't be tested
694 %! vopt = odeset ('OutputFcn', @fout, 'OutputSel', 1, 'Refine', 5);
695 %! vsol = ode45d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
696 %!test %# Stats must add further elements in vsol
697 %! vopt = odeset ('Stats', 'on');
698 %! vsol = ode45d (@fexp, [0 2.5], [1; 0], 1, [1; 0], vopt);
699 %! assert (isfield (vsol, 'stats'));
700 %! assert (isfield (vsol.stats, 'nsteps'));
701 %!test %# InitialStep option
702 %! vopt = odeset ('InitialStep', 1e-8);
703 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
704 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
705 %!test %# MaxStep option
706 %! vopt = odeset ('MaxStep', 1e-2);
707 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
708 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
709 %!test %# Events option add further elements in vsol
710 %! vopt = odeset ('Events', @feve);
711 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
712 %! assert (isfield (vsol, 'ie'));
713 %! assert (vsol.ie, [1; 1]);
714 %! assert (isfield (vsol, 'xe'));
715 %! assert (isfield (vsol, 'ye'));
716 %!test %# Events option, now stop integration
717 %! warning ('off', 'OdePkg:HideWarning');
718 %! vopt = odeset ('Events', @fevn, 'NormControl', 'on');
719 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
720 %! assert ([vsol.ie, vsol.xe, vsol.ye], ...
721 %! [1.0000, 2.9219, -0.2127, -0.2671], 0.5);
722 %!test %# Events option, five output arguments
723 %! vopt = odeset ('Events', @fevn, 'NormControl', 'on');
724 %! [vt, vy, vxe, vye, vie] = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
725 %! assert ([vie, vxe, vye], ...
726 %! [1.0000, 2.9219, -0.2127, -0.2671], 0.5);
728 %! %# test for Jacobian option is missing
729 %! %# test for Jacobian (being a sparse matrix) is missing
730 %! %# test for JPattern option is missing
731 %! %# test for Vectorized option is missing
732 %! %# test for NewtonTol option is missing
733 %! %# test for MaxNewtonIterations option is missing
735 %!test %# Mass option as function
736 %! vopt = odeset ('Mass', eye (2,2));
737 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
738 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
739 %!test %# Mass option as matrix
740 %! vopt = odeset ('Mass', eye (2,2));
741 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
742 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
743 %!test %# Mass option as sparse matrix
744 %! vopt = odeset ('Mass', sparse (eye (2,2)));
745 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
746 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
747 %!test %# Mass option as function and sparse matrix
748 %! vopt = odeset ('Mass', @fmsa);
749 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
750 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
751 %!test %# Mass option as function and MStateDependence
752 %! vopt = odeset ('Mass', @fmas, 'MStateDependence', 'strong');
753 %! vsol = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
754 %! assert ([vsol.x(end), vsol.y(end,:)], [5, fref], 0.5);
755 %!test %# Set BDF option to something else than default
756 %! vopt = odeset ('BDF', 'on');
757 %! [vt, vy] = ode45d (@fexp, [0 5], [1; 0], 1, [1; 0], vopt);
758 %! assert ([vt(end), vy(end,:)], [5, fref], 0.5);
760 %! %# test for MvPattern option is missing
761 %! %# test for InitialSlope option is missing
762 %! %# test for MaxOrder option is missing
764 %! warning ('on', 'OdePkg:InvalidOption');
766 %# Local Variables: ***