1 %# Copyright (C) 2007-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{solution}] =} odepkg_testsuite_chemakzo (@var{@@solver}, @var{reltol})
20 %# If this function is called with two input arguments and the first input argument @var{@@solver} is a function handle describing an OdePkg solver and the second input argument @var{reltol} is a double scalar describing the relative error tolerance then return a cell array @var{solution} with performance informations about the chemical AKZO Nobel testsuite of differential algebraic equations after solving (DAE--test).
22 %# Run examples with the command
24 %# demo odepkg_testsuite_chemakzo
27 %# This function has been ported from the "Test Set for IVP solvers" which is developed by the INdAM Bari unit project group "Codes and Test Problems for Differential Equations", coordinator F. Mazzia.
32 function vret = odepkg_testsuite_chemakzo (vhandle, vrtol)
34 if (nargin ~= 2) %# Check number and types of all input arguments
35 help ('odepkg_testsuite_chemakzo');
36 error ('OdePkg:InvalidArgument', ...
37 'Number of input arguments must be exactly two');
38 elseif (~isa (vhandle, 'function_handle') || ~isscalar (vrtol))
42 vret{1} = vhandle; %# The handle for the solver that is used
43 vret{2} = vrtol; %# The value for the realtive tolerance
44 vret{3} = vret{2}; %# The value for the absolute tolerance
45 vret{4} = vret{2}; %# The value for the first time step
46 %# Write a debug message on the screen, because this testsuite function
47 %# may be called more than once from a loop over all solvers present
48 fprintf (1, ['Testsuite AKZO, testing solver %7s with relative', ...
49 ' tolerance %2.0e\n'], func2str (vret{1}), vrtol); fflush (1);
51 %# Setting the integration algorithms option values
52 vstart = 0.0; %# The point of time when solving is started
53 vstop = 180.0; %# The point of time when solving is stoped
54 vinit = odepkg_testsuite_chemakzoinit; %# The initial values
55 vmass = odepkg_testsuite_chemakzomass; %# The mass matrix
57 vopt = odeset ('Refine', 0, 'RelTol', vret{2}, 'AbsTol', vret{3}, ...
58 'InitialStep', vret{4}, 'Stats', 'on', 'NormControl', 'off', ...
59 'Jacobian', @odepkg_testsuite_chemakzojac, 'Mass', vmass, ...
60 'MaxStep', vstop-vstart); %# , 'OutputFcn', @odeplot);
62 %# Calculate the algorithm, start timer and do solving
63 tic; vsol = feval (vhandle, @odepkg_testsuite_chemakzofun, ...
64 [vstart, vstop], vinit, vopt);
65 vret{12} = toc; %# The value for the elapsed time
66 vref = odepkg_testsuite_chemakzoref; %# Get the reference solution vector
67 if (max (size (vsol.y(end,:))) == max (size (vref))), vlst = vsol.y(end,:);
68 elseif (max (size (vsol.y(:,end))) == max (size (vref))), vlst = vsol.y(:,end);
70 vret{5} = odepkg_testsuite_calcmescd (vlst, vref, vret{3}, vret{2});
71 vret{6} = odepkg_testsuite_calcscd (vlst, vref, vret{3}, vret{2});
72 vret{7} = vsol.stats.nsteps + vsol.stats.nfailed; %# The value for all evals
73 vret{8} = vsol.stats.nsteps; %# The value for success evals
74 vret{9} = vsol.stats.nfevals; %# The value for fun calls
75 vret{10} = vsol.stats.npds; %# The value for partial derivations
76 vret{11} = vsol.stats.ndecomps; %# The value for LU decompositions
78 %# Returns the results for the for the chemical AKZO problem
79 function f = odepkg_testsuite_chemakzofun (t, y, varargin)
80 k1 = 18.7; k2 = 0.58; k3 = 0.09; k4 = 0.42;
81 kbig = 34.4; kla = 3.3; ks = 115.83; po2 = 0.9;
85 %# error ('odepkg_testsuite_chemakzojac: Second input argument is negative');
88 r1 = k1 * y(1)^4 * sqrt (y(2));
89 r2 = k2 * y(3) * y(4);
90 r3 = k2 / kbig * y(1) * y(5);
91 r4 = k3 * y(1) * y(4)^2;
92 r5 = k4 * y(6)^2 * sqrt (y(2));
93 fin = kla * (po2 / hen - y(2));
95 f(1,1) = -2 * r1 + r2 - r3 - r4;
96 f(2,1) = -0.5 * r1 - r4 - 0.5 * r5 + fin;
97 f(3,1) = r1 - r2 + r3;
98 f(4,1) = - r2 + r3 - 2 * r4;
99 f(5,1) = r2 - r3 + r5;
100 f(6,1) = ks * y(1) * y(4) - y(6);
102 %# Returns the INITIAL values for the chemical AKZO problem
103 function vinit = odepkg_testsuite_chemakzoinit ()
104 vinit = [0.444, 0.00123, 0, 0.007, 0, 115.83 * 0.444 * 0.007];
106 %# Returns the JACOBIAN matrix for the chemical AKZO problem
107 function dfdy = odepkg_testsuite_chemakzojac (t, y, varargin)
108 k1 = 18.7; k2 = 0.58; k3 = 0.09; k4 = 0.42;
109 kbig = 34.4; kla = 3.3; ks = 115.83; po2 = 0.9;
113 %# error ('odepkg_testsuite_chemakzojac: Second input argument is negative');
117 r11 = 4 * k1 * y(1)^3 * sqrt (y(2));
118 r12 = 0.5 * k1 * y(1)^4 / sqrt (y(2));
121 r31 = (k2 / kbig) * y(5);
122 r35 = (k2 / kbig) * y(1);
124 r44 = 2 * k3 * y(1) * y(4);
125 r52 = 0.5 * k4 * y(6)^2 / sqrt (y(2));
126 r56 = 2 * k4 * y(6) * sqrt (y(2));
129 dfdy(1,1) = -2 * r11 - r31 - r41;
130 dfdy(1,2) = -2 * r12;
132 dfdy(1,4) = r24 - r44;
134 dfdy(2,1) = -0.5 * r11 - r41;
135 dfdy(2,2) = -0.5 * r12 - 0.5 * r52 + fin2;
137 dfdy(2,6) = -0.5 * r56;
138 dfdy(3,1) = r11 + r31;
143 dfdy(4,1) = r31 - 2 * r41;
145 dfdy(4,4) = -r24 - 2 * r44;
153 dfdy(6,1) = ks * y(4);
154 dfdy(6,4) = ks * y(1);
157 %# Returns the MASS matrix for the chemical AKZO problem
158 function mass = odepkg_testsuite_chemakzomass (t, y, varargin)
159 mass = [ 1, 0, 0, 0, 0, 0;
166 %# Returns the REFERENCE values for the chemical AKZO problem
167 function y = odepkg_testsuite_chemakzoref ()
168 y(1,1) = 0.11507949206617e+0;
169 y(2,1) = 0.12038314715677e-2;
170 y(3,1) = 0.16115628874079e+0;
171 y(4,1) = 0.36561564212492e-3;
172 y(5,1) = 0.17080108852644e-1;
173 y(6,1) = 0.48735313103074e-2;
176 %! vsolver = {@odebda, @oders, @ode2r, @ode5r, @odesx};
177 %! for vcnt=1:length (vsolver)
178 %! vakzo{vcnt,1} = odepkg_testsuite_chemakzo (vsolver{vcnt}, 1e-7);
182 %# Local Variables: ***