1 %# Copyright (C) 2006-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{}] =} ode54 (@var{@@fun}, @var{slot}, @var{init}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
19 %# @deftypefnx {Command} {[@var{sol}] =} ode54 (@var{@@fun}, @var{slot}, @var{init}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
20 %# @deftypefnx {Command} {[@var{t}, @var{y}, [@var{xe}, @var{ye}, @var{ie}]] =} ode54 (@var{@@fun}, @var{slot}, @var{init}, [@var{opt}], [@var{par1}, @var{par2}, @dots{}])
22 %# This function file can be used to solve a set of non--stiff ordinary differential equations (non--stiff ODEs) or non--stiff differential algebraic equations (non--stiff DAEs) with the well known explicit Runge--Kutta method of order (5,4).
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 ODEs 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{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 %# If this function is called with one return argument then return the solution @var{sol} of type structure array after solving the set of ODEs. 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}.
28 %# 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.
30 %# For example, solve an anonymous implementation of the Van der Pol equation
33 %# fvdb = @@(vt,vy) [vy(2); (1 - vy(1)^2) * vy(2) - vy(1)];
35 %# vopt = odeset ("RelTol", 1e-3, "AbsTol", 1e-3, \
36 %# "NormControl", "on", "OutputFcn", @@odeplot);
37 %# ode54 (fvdb, [0 20], [2 0], vopt);
44 %# 20010703 the function file "ode54.m" was written by Marc Compere
45 %# under the GPL for the use with this software. This function has been
46 %# taken as a base for the following implementation.
47 %# 20060810, Thomas Treichl
48 %# This function was adapted to the new syntax that is used by the
49 %# new OdePkg for Octave. An equivalent function in Matlab does not
52 function [varargout] = ode54 (vfun, vslot, vinit, varargin)
54 if (nargin == 0) %# Check number and types of all input arguments
56 error ('OdePkg:InvalidArgument', ...
57 'Number of input arguments must be greater than zero');
62 elseif ~(isa (vfun, 'function_handle') || isa (vfun, 'inline'))
63 error ('OdePkg:InvalidArgument', ...
64 'First input argument must be a valid function handle');
66 elseif (~isvector (vslot) || length (vslot) < 2)
67 error ('OdePkg:InvalidArgument', ...
68 'Second input argument must be a valid vector');
70 elseif (~isvector (vinit) || ~isnumeric (vinit))
71 error ('OdePkg:InvalidArgument', ...
72 'Third input argument must be a valid numerical value');
76 if (~isstruct (varargin{1}))
77 %# varargin{1:len} are parameters for vfun
79 vfunarguments = varargin;
81 elseif (length (varargin) > 1)
82 %# varargin{1} is an OdePkg options structure vopt
83 vodeoptions = odepkg_structure_check (varargin{1}, 'ode54');
84 vfunarguments = {varargin{2:length(varargin)}};
86 else %# if (isstruct (varargin{1}))
87 vodeoptions = odepkg_structure_check (varargin{1}, 'ode54');
92 else %# if (nargin == 3)
97 %# Start preprocessing, have a look which options are set in
98 %# vodeoptions, check if an invalid or unused option is set
99 vslot = vslot(:).'; %# Create a row vector
100 vinit = vinit(:).'; %# Create a row vector
101 if (length (vslot) > 2) %# Step size checking
102 vstepsizefixed = true;
104 vstepsizefixed = false;
107 %# Get the default options that can be set with 'odeset' temporarily
110 %# Implementation of the option RelTol has been finished. This option
111 %# can be set by the user to another value than default value.
112 if (isempty (vodeoptions.RelTol) && ~vstepsizefixed)
113 vodeoptions.RelTol = 1e-6;
114 warning ('OdePkg:InvalidArgument', ...
115 'Option "RelTol" not set, new value %f is used', vodeoptions.RelTol);
116 elseif (~isempty (vodeoptions.RelTol) && vstepsizefixed)
117 warning ('OdePkg:InvalidArgument', ...
118 'Option "RelTol" will be ignored if fixed time stamps are given');
121 %# Implementation of the option AbsTol has been finished. This option
122 %# can be set by the user to another value than default value.
123 if (isempty (vodeoptions.AbsTol) && ~vstepsizefixed)
124 vodeoptions.AbsTol = 1e-6;
125 warning ('OdePkg:InvalidArgument', ...
126 'Option "AbsTol" not set, new value %f is used', vodeoptions.AbsTol);
127 elseif (~isempty (vodeoptions.AbsTol) && vstepsizefixed)
128 warning ('OdePkg:InvalidArgument', ...
129 'Option "AbsTol" will be ignored if fixed time stamps are given');
131 vodeoptions.AbsTol = vodeoptions.AbsTol(:); %# Create column vector
134 %# Implementation of the option NormControl has been finished. This
135 %# option can be set by the user to another value than default value.
136 if (strcmp (vodeoptions.NormControl, 'on')) vnormcontrol = true;
137 else vnormcontrol = false; end
139 %# Implementation of the option NonNegative has been finished. This
140 %# option can be set by the user to another value than default value.
141 if (~isempty (vodeoptions.NonNegative))
142 if (isempty (vodeoptions.Mass)), vhavenonnegative = true;
144 vhavenonnegative = false;
145 warning ('OdePkg:InvalidArgument', ...
146 'Option "NonNegative" will be ignored if mass matrix is set');
148 else vhavenonnegative = false;
151 %# Implementation of the option OutputFcn has been finished. This
152 %# option can be set by the user to another value than default value.
153 if (isempty (vodeoptions.OutputFcn) && nargout == 0)
154 vodeoptions.OutputFcn = @odeplot;
155 vhaveoutputfunction = true;
156 elseif (isempty (vodeoptions.OutputFcn)), vhaveoutputfunction = false;
157 else vhaveoutputfunction = true;
160 %# Implementation of the option OutputSel has been finished. This
161 %# option can be set by the user to another value than default value.
162 if (~isempty (vodeoptions.OutputSel)), vhaveoutputselection = true;
163 else vhaveoutputselection = false; end
165 %# Implementation of the option OutputSave has been finished. This
166 %# option can be set by the user to another value than default value.
167 if (isempty (vodeoptions.OutputSave)), vodeoptions.OutputSave = 1;
170 %# Implementation of the option Refine has been finished. This option
171 %# can be set by the user to another value than default value.
172 if (vodeoptions.Refine > 0), vhaverefine = true;
173 else vhaverefine = false; end
175 %# Implementation of the option Stats has been finished. This option
176 %# can be set by the user to another value than default value.
178 %# Implementation of the option InitialStep has been finished. This
179 %# option can be set by the user to another value than default value.
180 if (isempty (vodeoptions.InitialStep) && ~vstepsizefixed)
181 vodeoptions.InitialStep = (vslot(1,2) - vslot(1,1)) / 10;
182 vodeoptions.InitialStep = vodeoptions.InitialStep / 10^vodeoptions.Refine;
183 warning ('OdePkg:InvalidArgument', ...
184 'Option "InitialStep" not set, new value %f is used', vodeoptions.InitialStep);
187 %# Implementation of the option MaxStep has been finished. This option
188 %# can be set by the user to another value than default value.
189 if (isempty (vodeoptions.MaxStep) && ~vstepsizefixed)
190 vodeoptions.MaxStep = abs (vslot(1,2) - vslot(1,1)) / 10;
191 warning ('OdePkg:InvalidArgument', ...
192 'Option "MaxStep" not set, new value %f is used', vodeoptions.MaxStep);
195 %# Implementation of the option Events has been finished. This option
196 %# can be set by the user to another value than default value.
197 if (~isempty (vodeoptions.Events)), vhaveeventfunction = true;
198 else vhaveeventfunction = false; end
200 %# The options 'Jacobian', 'JPattern' and 'Vectorized' will be ignored
201 %# by this solver because this solver uses an explicit Runge-Kutta
202 %# method and therefore no Jacobian calculation is necessary
203 if (~isequal (vodeoptions.Jacobian, vodetemp.Jacobian))
204 warning ('OdePkg:InvalidArgument', ...
205 'Option "Jacobian" will be ignored by this solver');
207 if (~isequal (vodeoptions.JPattern, vodetemp.JPattern))
208 warning ('OdePkg:InvalidArgument', ...
209 'Option "JPattern" will be ignored by this solver');
211 if (~isequal (vodeoptions.Vectorized, vodetemp.Vectorized))
212 warning ('OdePkg:InvalidArgument', ...
213 'Option "Vectorized" will be ignored by this solver');
215 if (~isequal (vodeoptions.NewtonTol, vodetemp.NewtonTol))
216 warning ('OdePkg:InvalidArgument', ...
217 'Option "NewtonTol" will be ignored by this solver');
219 if (~isequal (vodeoptions.MaxNewtonIterations,...
220 vodetemp.MaxNewtonIterations))
221 warning ('OdePkg:InvalidArgument', ...
222 'Option "MaxNewtonIterations" will be ignored by this solver');
225 %# Implementation of the option Mass has been finished. This option
226 %# can be set by the user to another value than default value.
227 if (~isempty (vodeoptions.Mass) && isnumeric (vodeoptions.Mass))
228 vhavemasshandle = false; vmass = vodeoptions.Mass; %# constant mass
229 elseif (isa (vodeoptions.Mass, 'function_handle'))
230 vhavemasshandle = true; %# mass defined by a function handle
231 else %# no mass matrix - creating a diag-matrix of ones for mass
232 vhavemasshandle = false; %# vmass = diag (ones (length (vinit), 1), 0);
235 %# Implementation of the option MStateDependence has been finished.
236 %# This option can be set by the user to another value than default
238 if (strcmp (vodeoptions.MStateDependence, 'none'))
239 vmassdependence = false;
240 else vmassdependence = true;
243 %# Other options that are not used by this solver. Print a warning
244 %# message to tell the user that the option(s) is/are ignored.
245 if (~isequal (vodeoptions.MvPattern, vodetemp.MvPattern))
246 warning ('OdePkg:InvalidArgument', ...
247 'Option "MvPattern" will be ignored by this solver');
249 if (~isequal (vodeoptions.MassSingular, vodetemp.MassSingular))
250 warning ('OdePkg:InvalidArgument', ...
251 'Option "MassSingular" will be ignored by this solver');
253 if (~isequal (vodeoptions.InitialSlope, vodetemp.InitialSlope))
254 warning ('OdePkg:InvalidArgument', ...
255 'Option "InitialSlope" will be ignored by this solver');
257 if (~isequal (vodeoptions.MaxOrder, vodetemp.MaxOrder))
258 warning ('OdePkg:InvalidArgument', ...
259 'Option "MaxOrder" will be ignored by this solver');
261 if (~isequal (vodeoptions.BDF, vodetemp.BDF))
262 warning ('OdePkg:InvalidArgument', ...
263 'Option "BDF" will be ignored by this solver');
266 %# Starting the initialisation of the core solver ode54
267 vtimestamp = vslot(1,1); %# timestamp = start time
268 vtimelength = length (vslot); %# length needed if fixed steps
269 vtimestop = vslot(1,vtimelength); %# stop time = last value
270 %# 20110611, reported by Nils Strunk
271 %# Make it possible to solve equations from negativ to zero,
272 %# eg. vres = ode54 (@(t,y) y, [-2 0], 2);
273 vdirection = sign (vtimestop - vtimestamp); %# Direction flag
276 if (sign (vodeoptions.InitialStep) == vdirection)
277 vstepsize = vodeoptions.InitialStep;
278 else %# Fix wrong direction of InitialStep.
279 vstepsize = - vodeoptions.InitialStep;
281 vminstepsize = (vtimestop - vtimestamp) / (1/eps);
282 else %# If step size is given then use the fixed time steps
283 vstepsize = vslot(1,2) - vslot(1,1);
284 vminstepsize = sign (vstepsize) * eps;
287 vretvaltime = vtimestamp; %# first timestamp output
288 vretvalresult = vinit; %# first solution output
290 %# Initialize the OutputFcn
291 if (vhaveoutputfunction)
292 if (vhaveoutputselection) vretout = vretvalresult(vodeoptions.OutputSel);
293 else vretout = vretvalresult; end
294 feval (vodeoptions.OutputFcn, vslot.', ...
295 vretout.', 'init', vfunarguments{:});
298 %# Initialize the EventFcn
299 if (vhaveeventfunction)
300 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
301 vretvalresult.', 'init', vfunarguments{:});
304 vpow = 1/5; %# 20071016, reported by Luis Randez
305 va = [0, 0, 0, 0, 0, 0; %# The Dormand-Prince 5(4) coefficients
306 1/5, 0, 0, 0, 0, 0; %# Coefficients proved on 20060827
307 3/40, 9/40, 0, 0, 0, 0; %# See p.91 in Ascher & Petzold
308 44/45, -56/15, 32/9, 0, 0, 0;
309 19372/6561, -25360/2187, 64448/6561, -212/729, 0, 0;
310 9017/3168, -355/33, 46732/5247, 49/176, -5103/18656, 0;
311 35/384, 0, 500/1113, 125/192, -2187/6784, 11/84];
312 %# 4th and 5th order b-coefficients
313 vb4 = [5179/57600; 0; 7571/16695; 393/640; -92097/339200; 187/2100; 1/40];
314 vb5 = [35/384; 0; 500/1113; 125/192; -2187/6784; 11/84; 0];
317 %# The solver main loop - stop if the endpoint has been reached
318 vcntloop = 2; vcntcycles = 1; vu = vinit; vk = vu.' * zeros(1,7);
319 vcntiter = 0; vunhandledtermination = true; vcntsave = 2;
321 %# FSAL optimizations: sent by Bruce Minaker 20110326, find the slope k1
322 %# (k1 is copied from last k7, so set k7) for first time step only
323 if (vhavemasshandle) %# Handle only the dynamic mass matrix,
324 if (vmassdependence) %# constant mass matrices have already
325 vmass = feval ... %# been set before (if any)
326 (vodeoptions.Mass, vtimestamp, vinit.', vfunarguments{:});
327 else %# if (vmassdependence == false)
328 vmass = feval ... %# then we only have the time argument
329 (vodeoptions.Mass, vtimestamp, vfunarguments{:});
331 vk(:,7) = vmass \ feval ...
332 (vfun, vtimestamp, vinit.', vfunarguments{:});
335 (vfun, vtimestamp, vinit.', vfunarguments{:});
338 while ((vdirection * (vtimestamp) < vdirection * (vtimestop)) && ...
339 (vdirection * (vstepsize) >= vdirection * (vminstepsize)))
341 %# Hit the endpoint of the time slot exactely
342 if (vdirection * (vtimestamp + vstepsize) > vdirection * vtimestop)
343 %# vstepsize = vtimestop - vdirection * vtimestamp;
344 %# 20110611, reported by Nils Strunk
345 %# The endpoint of the time slot must be hit exactly,
346 %# eg. vsol = ode54 (@(t,y) y, [0 -1], 1);
347 vstepsize = vdirection * abs (abs (vtimestop) - abs (vtimestamp));
350 %# Estimate the seven results when using this solver (FSAL)
351 %# skip the first result as we already know it from last time step
353 for j = 2:7 %# Start at two instead of one (FSAL)
354 vthetime = vtimestamp + vc(j,1) * vstepsize;
355 vtheinput = vu.' + vstepsize * vk(:,1:j-1) * va(j,1:j-1).';
356 if (vhavemasshandle) %# Handle only the dynamic mass matrix,
357 if (vmassdependence) %# constant mass matrices have already
358 vmass = feval ... %# been set before (if any)
359 (vodeoptions.Mass, vthetime, vtheinput, vfunarguments{:});
360 else %# if (vmassdependence == false)
361 vmass = feval ... %# then we only have the time argument
362 (vodeoptions.Mass, vthetime, vfunarguments{:});
364 vk(:,j) = vmass \ feval ...
365 (vfun, vthetime, vtheinput, vfunarguments{:});
368 (vfun, vthetime, vtheinput, vfunarguments{:});
372 %# Compute the 4th and the 5th order estimation
373 y4 = vu.' + vstepsize * (vk * vb4);
375 %# y5 = vu.' + vstepsize * (vk * vb5); vb5 is the same as va(6,:),
376 %# means that we already know y5 from the first six vk's (FSAL)
378 if (vhavenonnegative)
379 vu(vodeoptions.NonNegative) = abs (vu(vodeoptions.NonNegative));
380 y4(vodeoptions.NonNegative) = abs (y4(vodeoptions.NonNegative));
381 y5(vodeoptions.NonNegative) = abs (y5(vodeoptions.NonNegative));
383 if (vhaveoutputfunction && vhaverefine)
384 vSaveVUForRefine = vu;
387 %# Calculate the absolute local truncation error and the acceptable error
390 vdelta = abs (y5 - y4);
391 vtau = max (vodeoptions.RelTol * abs (vu.'), vodeoptions.AbsTol);
393 vdelta = norm (y5 - y4, Inf);
394 vtau = max (vodeoptions.RelTol * max (norm (vu.', Inf), 1.0), ...
397 else %# if (vstepsizefixed == true)
398 vdelta = 1; vtau = 2;
401 %# If the error is acceptable then update the vretval variables
402 if (all (vdelta <= vtau))
403 vtimestamp = vtimestamp + vstepsize;
404 vu = y5.'; %# MC2001: the higher order estimation as "local extrapolation"
405 %# Save the solution every vodeoptions.OutputSave steps
406 if (mod (vcntloop-1,vodeoptions.OutputSave) == 0)
407 if (vhaveoutputselection)
408 vretvaltime(vcntsave,:) = vtimestamp;
409 vretvalresult(vcntsave,:) = vu(vodeoptions.OutputSel);
411 vretvaltime(vcntsave,:) = vtimestamp;
412 vretvalresult(vcntsave,:) = vu;
414 vcntsave = vcntsave + 1;
416 vcntloop = vcntloop + 1; vcntiter = 0;
418 %# Call plot only if a valid result has been found, therefore this
419 %# code fragment has moved here. Stop integration if plot function
421 if (vhaveoutputfunction)
422 for vcnt = 0:vodeoptions.Refine %# Approximation between told and t
423 if (vhaverefine) %# Do interpolation
424 vapproxtime = (vcnt + 1) * vstepsize / (vodeoptions.Refine + 2);
425 vapproxvals = vSaveVUForRefine.' + vapproxtime * (vk * vb5);
426 vapproxtime = (vtimestamp - vstepsize) + vapproxtime;
429 vapproxtime = vtimestamp;
431 if (vhaveoutputselection)
432 vapproxvals = vapproxvals(vodeoptions.OutputSel);
434 vpltret = feval (vodeoptions.OutputFcn, vapproxtime, ...
435 vapproxvals, [], vfunarguments{:});
436 if vpltret %# Leave refinement loop
440 if (vpltret) %# Leave main loop
441 vunhandledtermination = false;
446 %# Call event only if a valid result has been found, therefore this
447 %# code fragment has moved here. Stop integration if veventbreak is
449 if (vhaveeventfunction)
451 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
452 vu(:), [], vfunarguments{:});
453 if (~isempty (vevent{1}) && vevent{1} == 1)
454 vretvaltime(vcntloop-1,:) = vevent{3}(end,:);
455 vretvalresult(vcntloop-1,:) = vevent{4}(end,:);
456 vunhandledtermination = false; break;
461 %# If we're here, then we've overwritten the k7 that we
462 %# need to copy to k1 to redo the step Since we copy k7
463 %# into k1, we'll put the k1 back in k7 for now, until it's
464 %# copied back again (FSAL)
465 end %# If the error is acceptable ...
467 %# Update the step size for the next integration step
469 %# 2008-20120425, reported by Marco Caliari
470 %# vdelta cannot be negative (because of the absolute value that
471 %# has been introduced) but it could be 0, then replace the zeros
472 %# with the maximum value of vdelta
473 vdelta(find (vdelta == 0)) = max (vdelta);
474 %# It could happen that max (vdelta) == 0 (ie. that the original
475 %# vdelta was 0), in that case we double the previous vstepsize
476 vdelta(find (vdelta == 0)) = max (vtau) .* (0.4 ^ (1 / vpow));
479 vstepsize = min (vodeoptions.MaxStep, ...
480 min (0.8 * vstepsize * (vtau ./ vdelta) .^ vpow));
482 vstepsize = max (- vodeoptions.MaxStep, ...
483 max (0.8 * vstepsize * (vtau ./ vdelta) .^ vpow));
486 else %# if (vstepsizefixed)
487 if (vcntloop <= vtimelength)
488 vstepsize = vslot(vcntloop) - vslot(vcntloop-1);
489 else %# Get out of the main integration loop
494 %# Update counters that count the number of iteration cycles
495 vcntcycles = vcntcycles + 1; %# Needed for cost statistics
496 vcntiter = vcntiter + 1; %# Needed to find iteration problems
498 %# Stop solving because the last 1000 steps no successful valid
499 %# value has been found
500 if (vcntiter >= 5000)
501 error (['Solving has not been successful. The iterative', ...
502 ' integration loop exited at time t = %f before endpoint at', ...
503 ' tend = %f was reached. This happened because the iterative', ...
504 ' integration loop does not find a valid solution at this time', ...
505 ' stamp. Try to reduce the value of "InitialStep" and/or', ...
506 ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
511 %# Check if integration of the ode has been successful
512 if (vdirection * vtimestamp < vdirection * vtimestop)
513 if (vunhandledtermination == true)
514 error ('OdePkg:InvalidArgument', ...
515 ['Solving has not been successful. The iterative', ...
516 ' integration loop exited at time t = %f', ...
517 ' before endpoint at tend = %f was reached. This may', ...
518 ' happen if the stepsize grows smaller than defined in', ...
519 ' vminstepsize. Try to reduce the value of "InitialStep" and/or', ...
520 ' "MaxStep" with the command "odeset".\n'], vtimestamp, vtimestop);
522 warning ('OdePkg:InvalidArgument', ...
523 ['Solver has been stopped by a call of "break" in', ...
524 ' the main iteration loop at time t = %f before endpoint at', ...
525 ' tend = %f was reached. This may happen because the @odeplot', ...
526 ' function returned "true" or the @event function returned "true".'], ...
527 vtimestamp, vtimestop);
531 %# Postprocessing, do whatever when terminating integration algorithm
532 if (vhaveoutputfunction) %# Cleanup plotter
533 feval (vodeoptions.OutputFcn, vtimestamp, ...
534 vu.', 'done', vfunarguments{:});
536 if (vhaveeventfunction) %# Cleanup event function handling
537 odepkg_event_handle (vodeoptions.Events, vtimestamp, ...
538 vu.', 'done', vfunarguments{:});
540 %# Save the last step, if not already saved
541 if (mod (vcntloop-2,vodeoptions.OutputSave) ~= 0)
542 vretvaltime(vcntsave,:) = vtimestamp;
543 vretvalresult(vcntsave,:) = vu;
546 %# Print additional information if option Stats is set
547 if (strcmp (vodeoptions.Stats, 'on'))
549 vnsteps = vcntloop-2; %# vcntloop from 2..end
550 vnfailed = (vcntcycles-1)-(vcntloop-2)+1; %# vcntcycl from 1..end
551 vnfevals = 6*(vcntcycles-1)+1; %# 7 stages & one first step (FSAL)
552 vndecomps = 0; %# number of LU decompositions
553 vnpds = 0; %# number of partial derivatives
554 vnlinsols = 0; %# no. of solutions of linear systems
555 %# Print cost statistics if no output argument is given
557 vmsg = fprintf (1, 'Number of successful steps: %d\n', vnsteps);
558 vmsg = fprintf (1, 'Number of failed attempts: %d\n', vnfailed);
559 vmsg = fprintf (1, 'Number of function calls: %d\n', vnfevals);
565 if (nargout == 1) %# Sort output variables, depends on nargout
566 varargout{1}.x = vretvaltime; %# Time stamps are saved in field x
567 varargout{1}.y = vretvalresult; %# Results are saved in field y
568 varargout{1}.solver = 'ode54'; %# Solver name is saved in field solver
569 if (vhaveeventfunction)
570 varargout{1}.ie = vevent{2}; %# Index info which event occured
571 varargout{1}.xe = vevent{3}; %# Time info when an event occured
572 varargout{1}.ye = vevent{4}; %# Results when an event occured
575 varargout{1}.stats = struct;
576 varargout{1}.stats.nsteps = vnsteps;
577 varargout{1}.stats.nfailed = vnfailed;
578 varargout{1}.stats.nfevals = vnfevals;
579 varargout{1}.stats.npds = vnpds;
580 varargout{1}.stats.ndecomps = vndecomps;
581 varargout{1}.stats.nlinsols = vnlinsols;
583 elseif (nargout == 2)
584 varargout{1} = vretvaltime; %# Time stamps are first output argument
585 varargout{2} = vretvalresult; %# Results are second output argument
586 elseif (nargout == 5)
587 varargout{1} = vretvaltime; %# Same as (nargout == 2)
588 varargout{2} = vretvalresult; %# Same as (nargout == 2)
589 varargout{3} = []; %# LabMat doesn't accept lines like
590 varargout{4} = []; %# varargout{3} = varargout{4} = [];
592 if (vhaveeventfunction)
593 varargout{3} = vevent{3}; %# Time info when an event occured
594 varargout{4} = vevent{4}; %# Results when an event occured
595 varargout{5} = vevent{2}; %# Index info which event occured
600 %! # We are using the "Van der Pol" implementation for all tests that
601 %! # are done for this function. We also define a Jacobian, Events,
602 %! # pseudo-Mass implementation. For further tests we also define a
603 %! # reference solution (computed at high accuracy) and an OutputFcn
604 %!function [ydot] = fpol (vt, vy, varargin) %# The Van der Pol
605 %! ydot = [vy(2); (1 - vy(1)^2) * vy(2) - vy(1)];
606 %!function [vjac] = fjac (vt, vy, varargin) %# its Jacobian
607 %! vjac = [0, 1; -1 - 2 * vy(1) * vy(2), 1 - vy(1)^2];
608 %!function [vjac] = fjcc (vt, vy, varargin) %# sparse type
609 %! vjac = sparse ([0, 1; -1 - 2 * vy(1) * vy(2), 1 - vy(1)^2]);
610 %!function [vval, vtrm, vdir] = feve (vt, vy, varargin)
611 %! vval = fpol (vt, vy, varargin); %# We use the derivatives
612 %! vtrm = zeros (2,1); %# that's why component 2
613 %! vdir = ones (2,1); %# seems to not be exact
614 %!function [vval, vtrm, vdir] = fevn (vt, vy, varargin)
615 %! vval = fpol (vt, vy, varargin); %# We use the derivatives
616 %! vtrm = ones (2,1); %# that's why component 2
617 %! vdir = ones (2,1); %# seems to not be exact
618 %!function [vmas] = fmas (vt, vy)
619 %! vmas = [1, 0; 0, 1]; %# Dummy mass matrix for tests
620 %!function [vmas] = fmsa (vt, vy)
621 %! vmas = sparse ([1, 0; 0, 1]); %# A sparse dummy matrix
622 %!function [vref] = fref () %# The computed reference sol
623 %! vref = [0.32331666704577, -1.83297456798624];
624 %!function [vout] = fout (vt, vy, vflag, varargin)
625 %! if (regexp (char (vflag), 'init') == 1)
626 %! if (any (size (vt) ~= [2, 1])) error ('"fout" step "init"'); end
627 %! elseif (isempty (vflag))
628 %! if (any (size (vt) ~= [1, 1])) error ('"fout" step "calc"'); end
630 %! elseif (regexp (char (vflag), 'done') == 1)
631 %! if (any (size (vt) ~= [1, 1])) error ('"fout" step "done"'); end
632 %! else error ('"fout" invalid vflag');
635 %! %# Turn off output of warning messages for all tests, turn them on
636 %! %# again if the last test is called
637 %!error %# input argument number one
638 %! warning ('off', 'OdePkg:InvalidArgument');
639 %! B = ode54 (1, [0 25], [3 15 1]);
640 %!error %# input argument number two
641 %! B = ode54 (@fpol, 1, [3 15 1]);
642 %!error %# input argument number three
643 %! B = ode54 (@flor, [0 25], 1);
644 %!test %# one output argument
645 %! vsol = ode54 (@fpol, [0 2], [2 0]);
646 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
647 %! assert (isfield (vsol, 'solver'));
648 %! assert (vsol.solver, 'ode54');
649 %!test %# two output arguments
650 %! [vt, vy] = ode54 (@fpol, [0 2], [2 0]);
651 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-3);
652 %!test %# five output arguments and no Events
653 %! [vt, vy, vxe, vye, vie] = ode54 (@fpol, [0 2], [2 0]);
654 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-3);
655 %! assert ([vie, vxe, vye], []);
656 %!test %# anonymous function instead of real function
657 %! fvdb = @(vt,vy) [vy(2); (1 - vy(1)^2) * vy(2) - vy(1)];
658 %! vsol = ode54 (fvdb, [0 2], [2 0]);
659 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
660 %!test %# extra input arguments passed trhough
661 %! vsol = ode54 (@fpol, [0 2], [2 0], 12, 13, 'KL');
662 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
663 %!test %# empty OdePkg structure *but* extra input arguments
665 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt, 12, 13, 'KL');
666 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
667 %!error %# strange OdePkg structure
668 %! vopt = struct ('foo', 1);
669 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
670 %!test %# Solve vdp in fixed step sizes
671 %! vsol = ode54 (@fpol, [0:0.1:2], [2 0]);
672 %! assert (vsol.x(:), [0:0.1:2]');
673 %! assert (vsol.y(end,:), fref, 1e-3);
674 %!test %# Solve in backward direction starting at t=0
675 %! vref = [-1.205364552835178, 0.951542399860817];
676 %! vsol = ode54 (@fpol, [0 -2], [2 0]);
677 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-3);
678 %!test %# Solve in backward direction starting at t=2
679 %! vref = [-1.205364552835178, 0.951542399860817];
680 %! vsol = ode54 (@fpol, [2 -2], fref);
681 %! assert ([vsol.x(end), vsol.y(end,:)], [-2, vref], 1e-3);
682 %!test %# Solve another anonymous function in backward direction
683 %! vref = [-1, 0.367879437558975];
684 %! vsol = ode54 (@(t,y) y, [0 -1], 1);
685 %! assert ([vsol.x(end), vsol.y(end,:)], vref, 1e-3);
686 %!test %# Solve another anonymous function below zero
687 %! vref = [0, 14.77810590694212];
688 %! vsol = ode54 (@(t,y) y, [-2 0], 2);
689 %! assert ([vsol.x(end), vsol.y(end,:)], vref, 1e-3);
690 %!test %# AbsTol option
691 %! vopt = odeset ('AbsTol', 1e-5);
692 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
693 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
694 %!test %# AbsTol and RelTol option
695 %! vopt = odeset ('AbsTol', 1e-8, 'RelTol', 1e-8);
696 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
697 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
698 %!test %# RelTol and NormControl option -- higher accuracy
699 %! vopt = odeset ('RelTol', 1e-8, 'NormControl', 'on');
700 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
701 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-4);
702 %!test %# Keeps initial values while integrating
703 %! vopt = odeset ('NonNegative', 2);
704 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
705 %! assert ([vsol.x(end), vsol.y(end,:)], [2, 2, 0], 1e-1);
706 %!test %# Details of OutputSel and Refine can't be tested
707 %! vopt = odeset ('OutputFcn', @fout, 'OutputSel', 1, 'Refine', 5);
708 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
709 %!test %# Details of OutputSave can't be tested
710 %! vopt = odeset ('OutputSave', 1, 'OutputSel', 1);
711 %! vsla = ode54 (@fpol, [0 2], [2 0], vopt);
712 %! vopt = odeset ('OutputSave', 2);
713 %! vslb = ode54 (@fpol, [0 2], [2 0], vopt);
714 %! assert (length (vsla.x) > length (vslb.x))
715 %!test %# Stats must add further elements in vsol
716 %! vopt = odeset ('Stats', 'on');
717 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
718 %! assert (isfield (vsol, 'stats'));
719 %! assert (isfield (vsol.stats, 'nsteps'));
720 %!test %# InitialStep option
721 %! vopt = odeset ('InitialStep', 1e-8);
722 %! vsol = ode54 (@fpol, [0 0.2], [2 0], vopt);
723 %! assert ([vsol.x(2)-vsol.x(1)], [1e-8], 1e-9);
724 %!test %# MaxStep option
725 %! vopt = odeset ('MaxStep', 1e-2);
726 %! vsol = ode54 (@fpol, [0 0.2], [2 0], vopt);
727 %! assert ([vsol.x(5)-vsol.x(4)], [1e-2], 1e-3);
728 %!test %# Events option add further elements in vsol
729 %! vopt = odeset ('Events', @feve, 'InitialStep', 1e-6);
730 %! vsol = ode54 (@fpol, [0 10], [2 0], vopt);
731 %! assert (isfield (vsol, 'ie'));
732 %! assert (vsol.ie(1), 2);
733 %! assert (isfield (vsol, 'xe'));
734 %! assert (isfield (vsol, 'ye'));
735 %!test %# Events option, now stop integration
736 %! vopt = odeset ('Events', @fevn, 'InitialStep', 1e-6);
737 %! vsol = ode54 (@fpol, [0 10], [2 0], vopt);
738 %! assert ([vsol.ie, vsol.xe, vsol.ye], ...
739 %! [2.0, 2.496110, -0.830550, -2.677589], 1e-1);
740 %!test %# Events option, five output arguments
741 %! vopt = odeset ('Events', @fevn, 'InitialStep', 1e-6);
742 %! [vt, vy, vxe, vye, vie] = ode54 (@fpol, [0 10], [2 0], vopt);
743 %! assert ([vie, vxe, vye], ...
744 %! [2.0, 2.496110, -0.830550, -2.677589], 1e-1);
745 %!test %# Jacobian option
746 %! vopt = odeset ('Jacobian', @fjac);
747 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
748 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
749 %!test %# Jacobian option and sparse return value
750 %! vopt = odeset ('Jacobian', @fjcc);
751 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
752 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
754 %! %# test for JPattern option is missing
755 %! %# test for Vectorized option is missing
756 %! %# test for NewtonTol option is missing
757 %! %# test for MaxNewtonIterations option is missing
759 %!test %# Mass option as function
760 %! vopt = odeset ('Mass', @fmas);
761 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
762 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
763 %!test %# Mass option as matrix
764 %! vopt = odeset ('Mass', eye (2,2));
765 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
766 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
767 %!test %# Mass option as sparse matrix
768 %! vopt = odeset ('Mass', sparse (eye (2,2)));
769 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
770 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
771 %!test %# Mass option as function and sparse matrix
772 %! vopt = odeset ('Mass', @fmsa);
773 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
774 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
775 %!test %# Mass option as function and MStateDependence
776 %! vopt = odeset ('Mass', @fmas, 'MStateDependence', 'strong');
777 %! vsol = ode54 (@fpol, [0 2], [2 0], vopt);
778 %! assert ([vsol.x(end), vsol.y(end,:)], [2, fref], 1e-3);
779 %!test %# Set BDF option to something else than default
780 %! vopt = odeset ('BDF', 'on');
781 %! [vt, vy] = ode54 (@fpol, [0 2], [2 0], vopt);
782 %! assert ([vt(end), vy(end,:)], [2, fref], 1e-3);
784 %! %# test for MvPattern option is missing
785 %! %# test for InitialSlope option is missing
786 %! %# test for MaxOrder option is missing
788 %! warning ('on', 'OdePkg:InvalidArgument');
790 %# Local Variables: ***