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