1######################################################################## 2## 3## Copyright (C) 2006-2021 The Octave Project Developers 4## 5## See the file COPYRIGHT.md in the top-level directory of this 6## distribution or <https://octave.org/copyright/>. 7## 8## This file is part of Octave. 9## 10## Octave is free software: you can redistribute it and/or modify it 11## under the terms of the GNU General Public License as published by 12## the Free Software Foundation, either version 3 of the License, or 13## (at your option) any later version. 14## 15## Octave is distributed in the hope that it will be useful, but 16## WITHOUT ANY WARRANTY; without even the implied warranty of 17## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 18## GNU General Public License for more details. 19## 20## You should have received a copy of the GNU General Public License 21## along with Octave; see the file COPYING. If not, see 22## <https://www.gnu.org/licenses/>. 23## 24######################################################################## 25 26## -*- texinfo -*- 27## @deftypefn {} {[@var{t}, @var{y}] =} ode23s (@var{fun}, @var{trange}, @var{init}) 28## @deftypefnx {} {[@var{t}, @var{y}] =} ode23s (@var{fun}, @var{trange}, @var{init}, @var{ode_opt}) 29## @deftypefnx {} {[@var{t}, @var{y}] =} ode23s (@dots{}, @var{par1}, @var{par2}, @dots{}) 30## @deftypefnx {} {[@var{t}, @var{y}, @var{te}, @var{ye}, @var{ie}] =} ode23s (@dots{}) 31## @deftypefnx {} {@var{solution} =} ode23s (@dots{}) 32## 33## Solve a set of stiff Ordinary Differential Equations (stiff ODEs) with a 34## @nospell{Rosenbrock} method of order (2,3). 35## 36## @var{fun} is a function handle, inline function, or string containing the 37## name of the function that defines the ODE: @code{M y' = f(t,y)}. The 38## function must accept two inputs where the first is time @var{t} and the 39## second is a column vector of unknowns @var{y}. @var{M} is a constant mass 40## matrix, non-singular and possibly sparse. Set the field @qcode{"Mass"} in 41## @var{odeopts} using @var{odeset} to specify a mass matrix. 42## 43## @var{trange} specifies the time interval over which the ODE will be 44## evaluated. Typically, it is a two-element vector specifying the initial 45## and final times (@code{[tinit, tfinal]}). If there are more than two 46## elements then the solution will also be evaluated at these intermediate 47## time instances using an interpolation procedure of the same order as the 48## one of the solver. 49## 50## By default, @code{ode23s} uses an adaptive timestep with the 51## @code{integrate_adaptive} algorithm. The tolerance for the timestep 52## computation may be changed by using the options @qcode{"RelTol"} and 53## @qcode{"AbsTol"}. 54## 55## @var{init} contains the initial value for the unknowns. If it is a row 56## vector then the solution @var{y} will be a matrix in which each column is 57## the solution for the corresponding initial value in @var{init}. 58## 59## The optional fourth argument @var{ode_opt} specifies non-default options to 60## the ODE solver. It is a structure generated by @code{odeset}. 61## @code{ode23s} will ignore the following options: @qcode{"BDF"}, 62## @qcode{"InitialSlope"}, @qcode{"MassSingular"}, @qcode{"MStateDependence"}, 63## @qcode{"MvPattern"}, @qcode{"MaxOrder"}, @qcode{"Non-negative"}. 64## 65## The function typically returns two outputs. Variable @var{t} is a 66## column vector and contains the times where the solution was found. The 67## output @var{y} is a matrix in which each column refers to a different 68## unknown of the problem and each row corresponds to a time in @var{t}. If 69## @var{trange} specifies intermediate time steps, only those will be returned. 70## 71## The output can also be returned as a structure @var{solution} which has a 72## field @var{x} containing a row vector of times where the solution was 73## evaluated and a field @var{y} containing the solution matrix such that each 74## column corresponds to a time in @var{x}. Use 75## @w{@code{fieldnames (@var{solution})}} to see the other fields and 76## additional information returned. 77## 78## If using the @qcode{"Events"} option then three additional outputs may be 79## returned. @var{te} holds the time when an Event function returned a zero. 80## @var{ye} holds the value of the solution at time @var{te}. @var{ie} 81## contains an index indicating which Event function was triggered in the case 82## of multiple Event functions. 83## 84## Example: Solve the stiff @nospell{Van der Pol} equation 85## 86## @example 87## @group 88## f = @@(@var{t},@var{y}) [@var{y}(2); 1000*(1 - @var{y}(1)^2) * @var{y}(2) - @var{y}(1)]; 89## opt = odeset ('Mass', [1 0; 0 1], 'MaxStep', 1e-1); 90## [vt, vy] = ode23s (f, [0 2000], [2 0], opt); 91## @end group 92## @end example 93## @seealso{odeset, daspk, dassl} 94## @end deftypefn 95 96function varargout = ode23s (fun, trange, init, varargin) 97 98 if (nargin < 3) 99 print_usage (); 100 endif 101 102 solver = "ode23s"; 103 order = 2; 104 105 if (nargin >= 4) 106 if (! isstruct (varargin{1})) 107 ## varargin{1:len} are parameters for fun 108 odeopts = odeset (); 109 funarguments = varargin; 110 elseif (numel (varargin) > 1) 111 ## varargin{1} is an ODE options structure opt 112 odeopts = varargin{1}; 113 funarguments = {varargin{2:numel (varargin)}}; 114 else 115 ## varargin{1} is an ODE options structure opt 116 odeopts = varargin{1}; 117 funarguments = {}; 118 endif 119 else # nargin == 3 120 odeopts = odeset (); 121 funarguments = {}; 122 endif 123 124 if (! isnumeric (trange) || ! isvector (trange)) 125 error ("Octave:invalid-input-arg", 126 "ode23s: TRANGE must be a numeric vector"); 127 endif 128 129 if (numel (trange) < 2) 130 error ("Octave:invalid-input-arg", 131 "ode23s: TRANGE must contain at least 2 elements"); 132 elseif (trange(2) == trange(1)) 133 error ("Octave:invalid-input-arg", 134 "ode23s: invalid time span, TRANGE(1) == TRANGE(2)"); 135 else 136 direction = sign (trange(2) - trange(1)); 137 endif 138 trange = trange(:); 139 140 if (! isnumeric (init) || ! isvector (init)) 141 error ("Octave:invalid-input-arg", 142 "ode23s: INIT must be a numeric vector"); 143 endif 144 init = init(:); 145 146 if (ischar (fun)) 147 if (! exist (fun)) 148 error ("Octave:invalid-input-arg", 149 ['ode23s: function "' fun '" not found']); 150 endif 151 fun = str2func (fun); 152 endif 153 if (! is_function_handle (fun)) 154 error ("Octave:invalid-input-arg", 155 "ode23s: FUN must be a valid function handle"); 156 endif 157 158 ## Start preprocessing, have a look which options are set in odeopts, 159 ## check if an invalid or unused option is set. 160 [defaults, classes, attributes] = odedefaults (numel (init), 161 trange(1), trange(end)); 162 163 persistent ode23s_ignore_options = ... 164 {"BDF", "InitialSlope", "MassSingular", "MStateDependence", ... 165 "MvPattern", "MaxOrder", "NonNegative"}; 166 167 defaults = rmfield (defaults, ode23s_ignore_options); 168 classes = rmfield (classes, ode23s_ignore_options); 169 attributes = rmfield (attributes, ode23s_ignore_options); 170 171 odeopts = odemergeopts ("ode23s", odeopts, defaults, classes, attributes); 172 173 odeopts.funarguments = funarguments; 174 odeopts.direction = direction; 175 ## ode23s ignores "NonNegative" option, but integrate_adaptive needs it... 176 odeopts.havenonnegative = false; 177 178 if (isempty (odeopts.OutputFcn) && nargout == 0) 179 odeopts.OutputFcn = @odeplot; 180 odeopts.haveoutputfunction = true; 181 else 182 odeopts.haveoutputfunction = ! isempty (odeopts.OutputFcn); 183 endif 184 185 if (isempty (odeopts.InitialStep)) 186 odeopts.InitialStep = odeopts.direction * ... 187 starting_stepsize (order, fun, trange(1), init, 188 odeopts.AbsTol, odeopts.RelTol, 189 strcmpi (odeopts.NormControl, 190 "on"), 191 odeopts.funarguments); 192 endif 193 194 if (! isempty (odeopts.Mass) && isnumeric (odeopts.Mass)) 195 havemasshandle = false; 196 mass = odeopts.Mass; # constant mass 197 elseif (isa (odeopts.Mass, "function_handle")) 198 havemasshandle = true; # mass defined by a function handle 199 odeopts.Mass = feval (odeopts.Mass, trange(1), init, 200 odeopts.funarguments{:}); 201 else # no mass matrix - create a diag-matrix of ones for mass 202 havemasshandle = false; 203 odeopts.Mass = diag (ones (length (init), 1), 0); 204 endif 205 206 ## Starting the initialization of the core solver ode23s 207 208 if (nargout == 1) 209 ## Single output requires auto-selected intermediate times, 210 ## which is obtained by NOT specifying specific solution times. 211 trange = [trange(1); trange(end)]; 212 odeopts.Refine = []; # disable Refine when single output requested 213 elseif (numel (trange) > 2) 214 odeopts.Refine = []; # disable Refine when specific times requested 215 endif 216 217 solution = integrate_adaptive (@runge_kutta_23s, ... 218 order, fun, trange, init, odeopts); 219 220 ## Postprocessing, do whatever when terminating integration algorithm 221 if (odeopts.haveoutputfunction) # Cleanup plotter 222 feval (odeopts.OutputFcn, [], [], "done", odeopts.funarguments{:}); 223 endif 224 if (! isempty (odeopts.Events)) # Cleanup event function handling 225 ode_event_handler (odeopts.Events, solution.t(end), ... 226 solution.x(end,:).', "done", odeopts.funarguments{:}); 227 endif 228 229 ## Print additional information if option Stats is set 230 if (strcmpi (odeopts.Stats, "on")) 231 nsteps = solution.cntloop; # cntloop from 2..end 232 nfailed = solution.cntcycles - nsteps; # cntcycl from 1..end 233 nfevals = 5 * solution.cntcycles; # number of ode evaluations 234 ndecomps = nsteps; # number of LU decompositions 235 npds = 0; # number of partial derivatives 236 nlinsols = 3 * nsteps; # no. of solutions of linear systems 237 238 printf ("Number of successful steps: %d\n", nsteps); 239 printf ("Number of failed attempts: %d\n", nfailed); 240 printf ("Number of function calls: %d\n", nfevals); 241 endif 242 243 if (nargout == 2) 244 varargout{1} = solution.t; # Time stamps are first output argument 245 varargout{2} = solution.x; # Results are second output argument 246 elseif (nargout == 1) 247 varargout{1}.x = solution.t.'; # Time stamps saved in field x (row vector) 248 varargout{1}.y = solution.x.'; # Results are saved in field y (row vector) 249 varargout{1}.solver = solver; # Solver name is saved in field solver 250 if (! isempty (odeopts.Events)) 251 varargout{1}.xe = solution.event{3}; # Time info when an event occurred 252 varargout{1}.ye = solution.event{4}; # Results when an event occurred 253 varargout{1}.ie = solution.event{2}; # Index info which event occurred 254 endif 255 if (strcmpi (odeopts.Stats, "on")) 256 varargout{1}.stats = struct (); 257 varargout{1}.stats.nsteps = nsteps; 258 varargout{1}.stats.nfailed = nfailed; 259 varargout{1}.stats.nfevals = nfevals; 260 varargout{1}.stats.npds = npds; 261 varargout{1}.stats.ndecomps = ndecomps; 262 varargout{1}.stats.nlinsols = nlinsols; 263 endif 264 elseif (nargout == 5) 265 varargout = cell (1,5); 266 varargout{1} = solution.t; 267 varargout{2} = solution.x; 268 if (! isempty (odeopts.Events)) 269 varargout{3} = solution.event{3}; # Time info when an event occurred 270 varargout{4} = solution.event{4}; # Results when an event occurred 271 varargout{5} = solution.event{2}; # Index info which event occurred 272 endif 273 endif 274 275endfunction 276 277 278%!demo 279%! ## Demo function: stiff Van Der Pol equation 280%! fun = @(t,y) [y(2); 10*(1-y(1)^2)*y(2)-y(1)]; 281%! ## Calling ode23s method 282%! tic () 283%! [vt, vy] = ode23s (fun, [0 20], [2 0]); 284%! toc () 285%! ## Plotting the result 286%! plot(vt,vy(:,1),'-o'); 287 288%!demo 289%! ## Demo function: stiff Van Der Pol equation 290%! fun = @(t,y) [y(2); 10*(1-y(1)^2)*y(2)-y(1)]; 291%! ## Calling ode23s method 292%! odeopts = odeset ("Jacobian", @(t,y) [0 1; -20*y(1)*y(2)-1, 10*(1-y(1)^2)], 293%! "InitialStep", 1e-3) 294%! tic () 295%! [vt, vy] = ode23s (fun, [0 20], [2 0], odeopts); 296%! toc () 297%! ## Plotting the result 298%! plot(vt,vy(:,1),'-o'); 299 300%!demo 301%! ## Demo function: stiff Van Der Pol equation 302%! fun = @(t,y) [y(2); 100*(1-y(1)^2)*y(2)-y(1)]; 303%! ## Calling ode23s method 304%! odeopts = odeset ("InitialStep", 1e-4); 305%! tic () 306%! [vt, vy] = ode23s (fun, [0 200], [2 0]); 307%! toc () 308%! ## Plotting the result 309%! plot(vt,vy(:,1),'-o'); 310 311%!demo 312%! ## Demo function: stiff Van Der Pol equation 313%! fun = @(t,y) [y(2); 100*(1-y(1)^2)*y(2)-y(1)]; 314%! ## Calling ode23s method 315%! odeopts = odeset ("Jacobian", @(t,y) [0 1; -200*y(1)*y(2)-1, 100*(1-y(1)^2)], 316%! "InitialStep", 1e-4); 317%! tic () 318%! [vt, vy] = ode23s (fun, [0 200], [2 0], odeopts); 319%! toc () 320%! ## Plotting the result 321%! plot(vt,vy(:,1),'-o'); 322 323%!demo 324%! ## Demonstrate convergence order for ode23s 325%! tol = 1e-5 ./ 10.^[0:5]; 326%! for i = 1 : numel (tol) 327%! opt = odeset ("RelTol", tol(i), "AbsTol", realmin); 328%! [t, y] = ode23s (@(t, y) -y, [0, 1], 1, opt); 329%! h(i) = 1 / (numel (t) - 1); 330%! err(i) = norm (y .* exp (t) - 1, Inf); 331%! endfor 332%! 333%! ## Estimate order visually 334%! loglog (h, tol, "-ob", 335%! h, err, "-b", 336%! h, (h/h(end)) .^ 2 .* tol(end), "k--", 337%! h, (h/h(end)) .^ 3 .* tol(end), "k-"); 338%! axis tight 339%! xlabel ("h"); 340%! ylabel ("err(h)"); 341%! title ("Convergence plot for ode23s"); 342%! legend ("imposed tolerance", "ode23s (relative) error", 343%! "order 2", "order 3", "location", "northwest"); 344%! 345%! ## Estimate order numerically 346%! p = diff (log (err)) ./ diff (log (h)) 347 348%!test 349%! [vt, vy] = ode23s (@(t,y) t - y + 1, [0 10], [1]); 350%! assert ([vt(end), vy(end)], [10, exp(-10) + 10], 1e-3); 351 352%!test 353%! opts = odeset ('Mass', 5, 'Jacobian', -5, 'JConstant', 'on'); 354%! [vt, vy] = ode23s (@(t,y) 5 * (t - y + 1), [0 10], [1], opts); 355%! assert ([vt(end), vy(end)], [10, exp(-10) + 10], 1e-3); 356 357## We are using the "Van der Pol" implementation for all tests that are done 358## for this function. For further tests we also define a reference solution 359## (computed at high accuracy). 360%!function ydot = fpol (t, y) # The Van der Pol ODE 361%! ydot = [y(2); 10 * (1 - y(1)^2) * y(2) - y(1)]; 362%!endfunction 363%!function ydot = jac (t, y) # The Van der Pol ODE 364%! ydot = [0 1; -20 * y(1) * y(2) - 1, 10 * (1 - y(1)^2)]; 365%!endfunction 366%!function ref = fref () # The computed reference sol 367%! ref = [1.8610687248524305 -0.0753216319179125]; 368%!endfunction 369%!function [val, trm, dir] = feve (t, y, varargin) 370%! val = fpol (t, y, varargin); # We use the derivatives 371%! trm = zeros (2,1); # that's why component 2 372%! dir = ones (2,1); # does not seem to be exact 373%!endfunction 374%!function [val, trm, dir] = fevn (t, y, varargin) 375%! val = fpol (t, y, varargin); # We use the derivatives 376%! trm = ones (2,1); # that's why component 2 377%! dir = ones (2,1); # does not seem to be exact 378%!endfunction 379%!function mas = fmas (t, y, varargin) 380%! mas = [1, 0; 0, 1]; # Dummy mass matrix for tests 381%!endfunction 382%!function mas = fmsa (t, y, varargin) 383%! mas = sparse ([1, 0; 0, 1]); # A sparse dummy matrix 384%!endfunction 385%!function out = fout (t, y, flag, varargin) 386%! out = false; 387%! if (strcmp (flag, "init")) 388%! if (! isequal (size (t), [2, 1])) 389%! error ('fout: step "init"'); 390%! endif 391%! elseif (isempty (flag)) 392%! if (! isequal (size (t), [1, 1])) 393%! error ('fout: step "calc"'); 394%! endif 395%! elseif (strcmp (flag, "done")) 396%! if (! isempty (t)) 397%! warning ('fout: step "done"'); 398%! endif 399%! else 400%! error ("fout: invalid flag <%s>", flag); 401%! endif 402%!endfunction 403%! 404%!test # two output arguments 405%! [t, y] = ode23s (@fpol, [0 2], [2 0]); 406%! assert ([t(end), y(end,:)], [2, fref], 1e-3); 407%!test # anonymous function instead of real function 408%! fvdp = @(t,y) [y(2); 10 * (1 - y(1)^2) * y(2) - y(1)]; 409%! [t, y] = ode23s (fvdp, [0 2], [2 0]); 410%! assert ([t(end), y(end,:)], [2, fref], 1e-3); 411%!test # extra input arguments passed through 412%! [t, y] = ode23s (@fpol, [0 2], [2 0], 12, 13, "KL"); 413%! assert ([t(end), y(end,:)], [2, fref], 1e-3); 414%!test # empty OdePkg structure *but* extra input arguments 415%! opt = odeset; 416%! [t, y] = ode23s (@fpol, [0 2], [2 0], opt, 12, 13, "KL"); 417%! assert ([t(end), y(end,:)], [2, fref], 1e-2); 418%!test # InitialStep option 419%! opt = odeset ("InitialStep", 1e-8); 420%! [t, y] = ode23s (@fpol, [0 0.2], [2 0], opt); 421%! assert ([t(2)-t(1)], [1e-8], 1e-9); 422%!test # MaxStep option 423%! opt = odeset ("MaxStep", 1e-3); 424%! sol = ode23s (@fpol, [0 0.2], [2 0], opt); 425%! assert ([sol.x(5)-sol.x(4)], [1e-3], 1e-4); 426%!test # AbsTol option 427%! opt = odeset ("AbsTol", 1e-5); 428%! sol = ode23s (@fpol, [0 2], [2 0], opt); 429%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 430%!test # AbsTol and RelTol option 431%! opt = odeset ("AbsTol", 1e-8, "RelTol", 1e-8); 432%! sol = ode23s (@fpol, [0 2], [2 0], opt); 433%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 434%!test # RelTol and NormControl option -- higher accuracy 435%! opt = odeset ("RelTol", 1e-8, "NormControl", "on"); 436%! sol = ode23s (@fpol, [0 2], [2 0], opt); 437%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-4); 438%!test # Details of OutputSel and Refine can't be tested 439%! opt = odeset ("OutputFcn", @fout, "OutputSel", 1, "Refine", 5); 440%! sol = ode23s (@fpol, [0 2], [2 0], opt); 441%!test # Stats must add further elements in sol 442%! opt = odeset ("Stats", "on"); 443%! stat_str = evalc ("sol = ode23s (@fpol, [0 2], [2 0], opt);"); 444%! assert (strncmp (stat_str, "Number of successful steps:", 27)); 445%! assert (isfield (sol, "stats")); 446%! assert (isfield (sol.stats, "nsteps")); 447%!test # Events option add further elements in sol 448%! opt = odeset ("Events", @feve); 449%! sol = ode23s (@fpol, [0 10], [2 0], opt); 450%! assert (isfield (sol, "ie")); 451%! assert (sol.ie(1), 2); 452%! assert (isfield (sol, "xe")); 453%! assert (isfield (sol, "ye")); 454%!test # Mass option as function 455%! opt = odeset ("Mass", @fmas); 456%! sol = ode23s (@fpol, [0 2], [2 0], opt); 457%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 458%!test # Mass option as matrix 459%! opt = odeset ("Mass", eye (2,2)); 460%! sol = ode23s (@fpol, [0 2], [2 0], opt); 461%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 462%!test # Mass option as sparse matrix 463%! opt = odeset ("Mass", sparse (eye (2,2))); 464%! sol = ode23s (@fpol, [0 2], [2 0], opt); 465%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 466%!test # Mass option as function and sparse matrix 467%! opt = odeset ("Mass", @fmsa); 468%! sol = ode23s (@fpol, [0 2], [2 0], opt); 469%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 470%!test # Jacobian option as function 471%! opt = odeset ('Jacobian', @jac); 472%! sol = ode23s (@fpol, [0 2], [2 0], opt); 473%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 474%!testif HAVE_UMFPACK # Sparse Jacobian 475%! jac = @(t, y) sparse ([0 1; -20*y(1)*y(2)-1, 10*(1-y(1)^2)]); 476%! opt = odeset ('Jacobian', jac); 477%! sol = ode23s (@fpol, [0 2], [2 0], opt); 478%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 479%!testif HAVE_UMFPACK # Jpattern 480%! S = sparse ([0 1; 1 1]); 481%! opt = odeset ("Jpattern", S); 482%! sol = ode23s (@fpol, [0 2], [2 0], opt); 483%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3); 484 485## Note: The following options have no effect on this solver 486## therefore it makes no sense to test them here: 487## 488## "BDF" 489## "InitialSlope" 490## "MassSingular" 491## "MStateDependence" 492## "MaxOrder" 493## "MvPattern" 494## "NonNegative" 495 496%!test # Check that imaginary part of solution does not get inverted 497%! sol = ode23s (@(x,y) 1, [0 1], 1i); 498%! assert (imag (sol.y), ones (size (sol.y))) 499%! [x, y] = ode23s (@(x,y) 1, [0 1], 1i); 500%! assert (imag (y), ones (size (y))) 501 502## Test input validation 503%!error ode23s () 504%!error ode23s (1) 505%!error ode23s (1,2) 506%!error <TRANGE must be a numeric> ode23s (@fpol, {[0 25]}, [3 15 1]) 507%!error <TRANGE must be a .* vector> ode23s (@fpol, [0 25; 25 0], [3 15 1]) 508%!error <TRANGE must contain at least 2 elements> ode23s (@fpol, [1], [3 15 1]) 509%!error <invalid time span> ode23s (@fpol, [1 1], [3 15 1]) 510%!error <INIT must be a numeric> ode23s (@fpol, [0 25], {[3 15 1]}) 511%!error <INIT must be a .* vector> ode23s (@fpol, [0 25], [3 15 1; 3 15 1]) 512%!error <FUN must be a valid function handle> ode23s (1, [0 25], [3 15 1]) 513