File: ode45.m

package info (click to toggle)
octave 9.4.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 144,300 kB
  • sloc: cpp: 332,784; ansic: 77,239; fortran: 20,963; objc: 9,396; sh: 8,213; yacc: 4,925; lex: 4,389; perl: 1,544; java: 1,366; awk: 1,259; makefile: 648; xml: 189
file content (568 lines) | stat: -rw-r--r-- 22,787 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
########################################################################
##
## Copyright (C) 2006-2024 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {[@var{t}, @var{y}] =} ode45 (@var{fcn}, @var{trange}, @var{init})
## @deftypefnx {} {[@var{t}, @var{y}] =} ode45 (@var{fcn}, @var{trange}, @var{init}, @var{ode_opt})
## @deftypefnx {} {[@var{t}, @var{y}, @var{te}, @var{ye}, @var{ie}] =} ode45 (@dots{})
## @deftypefnx {} {@var{solution} =} ode45 (@dots{})
## @deftypefnx {} {} ode45 (@dots{})
##
## Solve a set of non-stiff Ordinary Differential Equations (non-stiff ODEs)
## with the well known explicit @nospell{Dormand-Prince} method of order 4.
##
## @var{fcn} is a function handle, inline function, or string containing the
## name of the function that defines the ODE: @code{y' = f(t,y)}.  The function
## must accept two inputs where the first is time @var{t} and the second is a
## column vector of unknowns @var{y}.
##
## @var{trange} specifies the time interval over which the ODE will be
## evaluated.  Typically, it is a two-element vector specifying the initial and
## final times (@code{[tinit, tfinal]}).  If there are more than two elements
## then the solution will also be evaluated at these intermediate time
## instances.
##
## By default, @code{ode45} uses an adaptive timestep with the
## @code{integrate_adaptive} algorithm.  The tolerance for the timestep
## computation may be changed by using the options @qcode{"RelTol"} and
## @qcode{"AbsTol"}.
##
## @var{init} contains the initial value for the unknowns.  If it is a row
## vector then the solution @var{y} will be a matrix in which each column is
## the solution for the corresponding initial value in @var{init}.
##
## The optional fourth argument @var{ode_opt} specifies non-default options to
## the ODE solver.  It is a structure generated by @code{odeset}.
##
## The function typically returns two outputs.  Variable @var{t} is a
## column vector and contains the times where the solution was found.  The
## output @var{y} is a matrix in which each column refers to a different
## unknown of the problem and each row corresponds to a time in @var{t}.
##
## The output can also be returned as a structure @var{solution} which has a
## field @var{x} containing a row vector of times where the solution was
## evaluated and a field @var{y} containing the solution matrix such that each
## column corresponds to a time in @var{x}.  Use
## @w{@code{fieldnames (@var{solution})}}@ to see the other fields and
## additional information returned.
##
## If no output arguments are requested, and no @qcode{"OutputFcn"} is
## specified in @var{ode_opt}, then the @qcode{"OutputFcn"} is set to
## @code{odeplot} and the results of the solver are plotted immediately.
##
## If using the @qcode{"Events"} option then three additional outputs may be
## returned.  @var{te} holds the time when an Event function returned a zero.
## @var{ye} holds the value of the solution at time @var{te}.  @var{ie}
## contains an index indicating which Event function was triggered in the case
## of multiple Event functions.
##
## Example: Solve the @nospell{Van der Pol} equation
##
## @example
## @group
## fvdp = @@(@var{t},@var{y}) [@var{y}(2); (1 - @var{y}(1)^2) * @var{y}(2) - @var{y}(1)];
## [@var{t},@var{y}] = ode45 (fvdp, [0, 20], [2, 0]);
## @end group
## @end example
## @seealso{odeset, odeget, ode23, ode15s}
## @end deftypefn

function varargout = ode45 (fcn, trange, init, varargin)

  if (nargin < 3)
    print_usage ();
  endif

  solver = "ode45";
  order  = 5;  # runge_kutta_45_dorpri uses local extrapolation

  if (nargin >= 4)
    if (! isstruct (varargin{1}))
      ## varargin{1:len} are parameters for fcn
      odeopts = odeset ();
      funarguments = varargin;
    elseif (numel (varargin) > 1)
      ## varargin{1} is an ODE options structure opt
      odeopts = varargin{1};
      funarguments = {varargin{2:numel (varargin)}};
    else
      ## varargin{1} is an ODE options structure opt
      odeopts = varargin{1};
      funarguments = {};
    endif
  else  # nargin == 3
    odeopts = odeset ();
    funarguments = {};
  endif

  if (! isnumeric (trange) || ! isvector (trange))
    error ("Octave:invalid-input-arg",
           "ode45: TRANGE must be a numeric vector");
  endif

  if (numel (trange) < 2)
    error ("Octave:invalid-input-arg",
           "ode45: TRANGE must contain at least 2 elements");
  elseif (trange(1) == trange(2))
    error ("Octave:invalid-input-arg",
           "ode45: invalid time span, TRANGE(1) == TRANGE(2)");
  else
    direction = sign (trange(2) - trange(1));
  endif
  trange = trange(:);

  if (! isnumeric (init) || ! isvector (init))
    error ("Octave:invalid-input-arg",
           "ode45: INIT must be a numeric vector");
  endif
  init = init(:);

  if (ischar (fcn))
    if (! exist (fcn))
      error ("Octave:invalid-input-arg",
             ['ode45: function "' fcn '" not found']);
    endif
    fcn = str2func (fcn);
  endif
  if (! is_function_handle (fcn))
    error ("Octave:invalid-input-arg",
           "ode45: FCN must be a valid function handle");
  endif

  ## FIXME: Warn user if ! isempty (funarguments)
  ## Not a documented behavior and may be deprecated


  ## Start preprocessing, have a look which options are set in odeopts,
  ## check if an invalid or unused option is set
  [defaults, classes, attributes] = odedefaults (numel (init),
                                                 trange(1), trange(end));

  defaults = odeset (defaults, "Refine", 4);

  persistent ode45_ignore_options = ...
    {"BDF", "InitialSlope", "Jacobian", "JPattern",
     "MassSingular", "MaxOrder", "MvPattern", "Vectorized"};

  defaults   = rmfield (defaults, ode45_ignore_options);
  classes    = rmfield (classes, ode45_ignore_options);
  attributes = rmfield (attributes, ode45_ignore_options);

  odeopts = odemergeopts ("ode45", odeopts, defaults, classes, attributes);

  odeopts.funarguments = funarguments;
  odeopts.direction    = direction;

  if (! isempty (odeopts.NonNegative))
    if (isempty (odeopts.Mass))
      odeopts.havenonnegative = true;
    else
      odeopts.havenonnegative = false;
      warning ("Octave:invalid-input-arg",
               ['ode45: option "NonNegative" is ignored', ...
                " when mass matrix is set\n"]);
    endif
  else
    odeopts.havenonnegative = false;
  endif

  if (isempty (odeopts.OutputFcn) && nargout == 0)
    odeopts.OutputFcn = @odeplot;
    odeopts.haveoutputfunction = true;
  else
    odeopts.haveoutputfunction = ! isempty (odeopts.OutputFcn);
  endif

  if (isempty (odeopts.InitialStep))
    odeopts.InitialStep = odeopts.direction * ...
                          starting_stepsize (order, fcn, trange(1), init,
                                             odeopts.AbsTol, odeopts.RelTol,
                                             strcmpi (odeopts.NormControl,
                                             "on"), odeopts.funarguments);
  endif

  if (! isempty (odeopts.Mass))
    if (isnumeric (odeopts.Mass))
      havemasshandle = false;
      mass = odeopts.Mass;  # constant mass
    elseif (is_function_handle (odeopts.Mass))
      havemasshandle = true;    # mass defined by a function handle
    else
      error ("Octave:invalid-input-arg",
             'ode45: "Mass" field must be a function handle or square matrix');
    endif
  else  # no mass matrix - create a diag-matrix of ones for mass
    havemasshandle = false;   # mass = diag (ones (length (init), 1), 0);
  endif

  ## Starting the initialization of the core solver ode45

  if (havemasshandle)   # Handle only the dynamic mass matrix,
    if (! strcmp (odeopts.MStateDependence, "none"))
      ### constant mass matrices have already
      mass = @(t,x) odeopts.Mass (t, x, odeopts.funarguments{:});
      fcn = @(t,x) mass (t, x, odeopts.funarguments{:}) ...
                   \ fcn (t, x, odeopts.funarguments{:});
    else
      mass = @(t) odeopts.Mass (t, odeopts.funarguments{:});
      fcn = @(t,x) mass (t, odeopts.funarguments{:}) ...
                   \ fcn (t, x, odeopts.funarguments{:});
    endif
  endif

  if (numel (trange) > 2)
    odeopts.Refine = [];  # disable Refine when specific times requested
  endif

  solution = integrate_adaptive (@runge_kutta_45_dorpri,
                                 order, fcn, trange, init, odeopts);

  ## Postprocessing, do whatever when terminating integration algorithm
  if (odeopts.haveoutputfunction)  # Cleanup plotter
    feval (odeopts.OutputFcn, [], [], "done", odeopts.funarguments{:});
  endif
  if (! isempty (odeopts.Events))   # Cleanup event function handling
    ode_event_handler ([], [], [], [], [], "done");
  endif

  ## Print additional information if option Stats is set
  if (strcmpi (odeopts.Stats, "on"))
    nsteps    = solution.cntloop;             # cntloop from 2..end
    nfailed   = solution.cntcycles - nsteps;  # cntcycl from 1..end
    nfevals   = 6 * solution.cntcycles + 1;   # number of ode evaluations
    ndecomps  = 0;  # number of LU decompositions
    npds      = 0;  # number of partial derivatives
    nlinsols  = 0;  # no. of linear systems solutions

    printf ("Number of successful steps: %d\n", nsteps);
    printf ("Number of failed attempts:  %d\n", nfailed);
    printf ("Number of function calls:   %d\n", nfevals);
  endif

  if (nargout == 2)
    varargout{1} = solution.output_t; # Time stamps are first output argument
    varargout{2} = solution.output_x; # Results are second output argument
  elseif (nargout == 1)
    varargout{1}.x = solution.ode_t.'; #Time stamps saved in field x (row vect.)
    varargout{1}.y = solution.ode_x.'; #Results are saved in field y (row vect.)
    varargout{1}.solver = solver; # Solver name is saved in field solver
    if (! isempty (odeopts.Events))
      varargout{1}.xe = solution.event{3}.'; # Time info when an event occurred
      varargout{1}.ye = solution.event{4}.'; # Results when an event occurred
      varargout{1}.ie = solution.event{2}.'; # Index info which event occurred
    endif
    if (strcmpi (odeopts.Stats, "on"))
      varargout{1}.stats = struct ();
      varargout{1}.stats.nsteps   = nsteps;
      varargout{1}.stats.nfailed  = nfailed;
      varargout{1}.stats.nfevals  = nfevals;
      varargout{1}.stats.npds     = npds;
      varargout{1}.stats.ndecomps = ndecomps;
      varargout{1}.stats.nlinsols = nlinsols;
    endif
  elseif (nargout > 2)
    varargout = cell (1,5);
    varargout{1} = solution.output_t;
    varargout{2} = solution.output_x;
    if (! isempty (odeopts.Events))
      varargout{3} = solution.event{3};  # Time info when an event occurred
      varargout{4} = solution.event{4};  # Results when an event occurred
      varargout{5} = solution.event{2};  # Index info which event occurred
    endif
  endif

endfunction


%!demo
%! ## Demonstrate convergence order for ode45
%! tol = 1e-5 ./ 10.^[0:8];
%! for i = 1 : numel (tol)
%!   opt = odeset ("RelTol", tol(i), "AbsTol", realmin);
%!   [t, y] = ode45 (@(t, y) -y, [0, 1], 1, opt);
%!   h(i) = 1 / (numel (t) - 1);
%!   err(i) = norm (y .* exp (t) - 1, Inf);
%! endfor
%!
%! ## Estimate order visually
%! loglog (h, tol, "-ob",
%!         h, err, "-b",
%!         h, (h/h(end)) .^ 4 .* tol(end), "k--",
%!         h, (h/h(end)) .^ 5 .* tol(end), "k-");
%! axis tight
%! xlabel ("h");
%! ylabel ("err(h)");
%! title ("Convergence plot for ode45");
%! legend ("imposed tolerance", "ode45 (relative) error",
%!         "order 4", "order 5", "location", "northwest");
%!
%! ## Estimate order numerically
%! p = diff (log (err)) ./ diff (log (h))

## We are using the Van der Pol equation for all tests.
## Further tests also define a reference solution (computed at high accuracy)
%!function ydot = fpol (t, y, varargin)  # The Van der Pol ODE
%!  ydot = [y(2); (1 - y(1)^2) * y(2) - y(1)];
%!endfunction
%!function ref = fref ()       # The computed reference solution
%!  ref = [0.32331666704577, -1.83297456798624];
%!endfunction
%!function [val, trm, dir] = feve (t, y, varargin)
%!  val = fpol (t, y, varargin{:});  # We use the derivatives
%!  trm = zeros (2,1);            # that's why component 2
%!  dir = ones (2,1);             # does not seem to be exact
%!endfunction
%!function [val, trm, dir] = fevn (t, y, varargin)
%!  val = fpol (t, y, varargin{:});  # We use the derivatives
%!  trm = ones (2,1);             # that's why component 2
%!  dir = ones (2,1);             # does not seem to be exact
%!endfunction
%!function mas = fmas (t, y, varargin)
%!  mas = [1, 0; 0, 1];           # Dummy mass matrix for tests
%!endfunction
%!function mas = fmsa (t, y, varargin)
%!  mas = sparse ([1, 0; 0, 1]);  # A sparse dummy matrix
%!endfunction
%!function out = fout (t, y, flag, varargin)
%!  out = false;
%!  if (strcmp (flag, "init"))
%!    if (! isequal (size (t), [2, 1]))
%!      error ('fout: step "init"');
%!    endif
%!  elseif (isempty (flag))
%!  # Multiple steps can be sent in one function call
%!    if (! isequal ( size (t), size (y)))
%!      error ('fout: step "calc"');
%!    endif
%!  elseif (strcmp (flag, "done"))
%!    if (! isempty (t))
%!      warning ('fout: step "done"');
%!    endif
%!  else
%!    error ("fout: invalid flag <%s>", flag);
%!  endif
%!endfunction
%!function stop_solve = OutputSel_test (t, y, flag, x)
%!  ## x == 1: select y(1)
%!  ## x == 2: select y(2)
%!  ## x == 3: select y([1,2])
%!  persistent y_last
%!  if (strcmp (flag, "init"))
%!    y_last = y;
%!    if (x == 1 || x == 2)
%!      assert (length (y) == 1);
%!    elseif (x == 3)
%!      assert (length (y) == 2);
%!    endif
%!  elseif (strcmp (flag, "done"))
%!    y_exp = fref ().';
%!    if (x < 3)
%!      assert (y_last, y_exp(x), 1e-4);
%!    else
%!      assert (y_last, y_exp, 1e-4);
%!    endif
%!  else  # flag == ""
%!    y_last = y(:,end);
%!    if (x == 1 || x == 2)
%!      assert (length (t) == length (y));
%!    else
%!      assert (2 * length (t) == length (y(:)));
%!    endif
%!  endif
%!  stop_solve = 0;
%!endfunction
%!
%!test  # two output arguments
%! [t, y] = ode45 (@fpol, [0 2], [2 0]);
%! assert ([t(end), y(end,:)], [2, fref], 1e-2);
%!test  # not too many steps
%! [t, y] = ode45 (@fpol, [0 2], [2 0], odeset("Refine",1));
%! assert (size (t) < 20);
%!test  # correct number of steps with Refine
%! [t1, y1] = ode45 (@fpol, [0 2], [2 0], odeset ("Refine", 1));
%! [t2, y2] = ode45 (@fpol, [0 2], [2 0], odeset ("Refine", 4));
%! [t3, y3] = ode45 (@fpol, [0 2], [2 0]); #default Refine=4
%! s = ode45 (@fpol, [0 2], [2 0], odeset ("Refine", 4));
%! assert (length (t2) == length (t3));
%! assert (length (t2) == 4*length (t1) - 3);
%! assert (length (s.x) == length (t1));
%!test  # anonymous function instead of real function
%! fvdp = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)];
%! [t, y] = ode45 (fvdp, [0 2], [2 0]);
%! assert ([t(end), y(end,:)], [2, fref], 1e-2);
%!test  # string instead of function
%! [t, y] = ode45 ("fpol", [0 2], [2 0]);
%! assert ([t(end), y(end,:)], [2, fref], 1e-2);
%!test  # extra input arguments passed through
%! [t, y] = ode45 (@fpol, [0 2], [2 0], 12, 13, "KL");
%! assert ([t(end), y(end,:)], [2, fref], 1e-2);
%!test  # empty ODEOPT structure *but* extra input arguments
%! opt = odeset;
%! [t, y] = ode45 (@fpol, [0 2], [2 0], opt, 12, 13, "KL");
%! assert ([t(end), y(end,:)], [2, fref], 1e-2);
%!test  # Solve another anonymous function below zero
%! vref = [0, 14.77810590694212];
%! [t, y] = ode45 (@(t,y) y, [-2 0], 2);
%! assert ([t(end), y(end,:)], vref, 1e-1);
%!test  # InitialStep option
%! opt = odeset ("InitialStep", 1e-8, "Refine", 1);
%! [t, y] = ode45 (@fpol, [0 0.2], [2 0], opt);
%! assert ([t(2)-t(1)], [1e-8], 1e-9);
%!test  # MaxStep option
%! opt = odeset ("MaxStep", 1e-3);
%! sol = ode45 (@fpol, [0 0.2], [2 0], opt);
%! assert ([sol.x(5)-sol.x(4)], [1e-3], 1e-3);
%!test  # Solve with intermediate step
%! [t, y] = ode45 (@fpol, [0 1 2], [2 0]);
%! assert (any ((t-1) == 0));
%! assert ([t(end), y(end,:)], [2, fref], 1e-3);
%!test  # Solve in backward direction starting at t=0
%! vref = [-1.205364552835178, 0.951542399860817];
%! sol = ode45 (@fpol, [0 -2], [2 0]);
%! assert ([sol.x(end); sol.y(:,end)], [-2; vref'], 1e-2);
%!test  # Solve in backward direction starting at t=2
%! vref = [-1.205364552835178, 0.951542399860817];
%! sol = ode45 (@fpol, [2 -2], fref);
%! assert ([sol.x(end); sol.y(:,end)], [-2; vref'], 1e-2);
%!test  # Solve in backward direction starting at t=2, with intermediate step
%! vref = [-1.205364552835178, 0.951542399860817];
%! [t, y] = ode45 (@fpol, [2 0 -2], fref);
%! idx = find (y < 0, 1, "first") - 1;
%! assert ([t(idx), y(idx,:)], [0,2,0], 1e-2);
%! assert ([t(end), y(end,:)], [-2, vref], 1e-2);
%!test  # Solve another anonymous function in backward direction
%! vref = [-1, 0.367879437558975];
%! sol = ode45 (@(t,y) y, [0 -1], 1);
%! assert ([sol.x(end); sol.y(:,end)], vref', 1e-3);
%!test  # Solve another anonymous function below zero
%! vref = [0, 14.77810590694212];
%! sol = ode45 (@(t,y) y, [-2 0], 2);
%! assert ([sol.x(end); sol.y(:,end)], vref', 1e-3);
%!test  # Solve in backward direction starting at t=0 with MaxStep option
%! vref = [-1.205364552835178, 0.951542399860817];
%! opt = odeset ("MaxStep", 1e-3);
%! sol = ode45 (@fpol, [0 -2], [2 0], opt);
%! assert ([abs(sol.x(8)-sol.x(7))], [1e-3], 1e-3);
%! assert ([sol.x(end); sol.y(:,end)], [-2; vref'], 1e-3);
%!test  # AbsTol option
%! opt = odeset ("AbsTol", 1e-5);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # AbsTol and RelTol option
%! opt = odeset ("AbsTol", 1e-8, "RelTol", 1e-8);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # RelTol and NormControl option -- higher accuracy
%! opt = odeset ("RelTol", 1e-8, "NormControl", "on");
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-5);
%!test  # Keeps initial values while integrating
%! opt = odeset ("NonNegative", 2);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; 2; 0], 0.5);
%!test  # OutputSel 1 (see function OutputSel_test for asserts)
%! opt = odeset ("OutputFcn", @(t, y, flag) OutputSel_test (t, y, flag, 1), ...
%!               "OutputSel", 1);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%!test  # OutputSel 2 (see function OutputSel_test for asserts)
%! opt = odeset ("OutputFcn", @(t, y, flag) OutputSel_test (t, y, flag, 2), ...
%!               "OutputSel", 2);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%!test  # OutputSel [1,2] (see function OutputSel_test for asserts)
%! opt = odeset ("OutputFcn", @(t, y, flag) OutputSel_test (t, y, flag, 3), ...
%!               "OutputSel", [1,2]);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%!test  # Stats must add further elements in sol
%! opt = odeset ("Stats", "on");
%! stat_str = evalc ("sol = ode45 (@fpol, [0 2], [2 0], opt);");
%! assert (strncmp (stat_str, "Number of successful steps:", 27));
%! assert (isfield (sol, "stats"));
%! assert (isfield (sol.stats, "nsteps"));
%!test  # Events option add further elements in sol
%! opt = odeset ("Events", @feve);
%! sol = ode45 (@fpol, [0 10], [2 0], opt);
%! assert (isfield (sol, "ie"));
%! assert (sol.ie(1), 2);
%! assert (isfield (sol, "xe"));
%! assert (isfield (sol, "ye"));
%!test  # Events option, now stop integration
%! opt = odeset ("Events", @fevn, "NormControl", "on");
%! sol = ode45 (@fpol, [0 10], [2 0], opt);
%! assert ([sol.ie, sol.xe, sol.ye.'],
%!         [2.0, 2.496110, -0.830550, -2.677589], 2e-3);
%!test  # Events option, five output arguments
%! opt = odeset ("Events", @fevn, "NormControl", "on");
%! [t, y, vxe, ye, vie] = ode45 (@fpol, [0 10], [2 0], opt);
%! assert ([vie, vxe, ye],
%!         [2.0, 2.496110, -0.830550, -2.677589], 2e-3);
%!test  # Mass option as function
%! opt = odeset ("Mass", @fmas);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # Mass option as matrix
%! opt = odeset ("Mass", eye (2,2));
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # Mass option as sparse matrix
%! opt = odeset ("Mass", sparse (eye (2,2)));
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # Mass option as function and sparse matrix
%! opt = odeset ("Mass", @fmsa);
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);
%!test  # Mass option as function and MStateDependence
%! opt = odeset ("Mass", @fmas, "MStateDependence", "strong");
%! sol = ode45 (@fpol, [0 2], [2 0], opt);
%! assert ([sol.x(end); sol.y(:,end)], [2; fref'], 1e-3);

## Note: The following options have no effect on this solver
##       therefore it makes no sense to test them here:
##
## "BDF"
## "InitialSlope"
## "JPattern"
## "Jacobian"
## "MassSingular"
## "MaxOrder"
## "MvPattern"
## "Vectorized"

%!test  # Check that imaginary part of solution does not get inverted
%! sol = ode45 (@(x,y) 1, [0 1], 1i);
%! assert (imag (sol.y), ones (size (sol.y)));
%! [x, y] = ode45 (@(x,y) 1, [0 1], 1i, odeset ("Refine", 1));
%! assert (imag (y), ones (size (y)));

%!error <Invalid call> ode45 ()
%!error <Invalid call> ode45 (1)
%!error <Invalid call> ode45 (1,2)
%!error <TRANGE must be a numeric> ode45 (@fpol, {[0 25]}, [3 15 1])
%!error <TRANGE must be a .* vector> ode45 (@fpol, [0 25; 25 0], [3 15 1])
%!error <TRANGE must contain at least 2 elements> ode45 (@fpol, [1], [3 15 1])
%!error <invalid time span> ode45 (@fpol, [1 1], [3 15 1])
%!error <INIT must be a numeric> ode45 (@fpol, [0 25], {[3 15 1]})
%!error <INIT must be a .* vector> ode45 (@fpol, [0 25], [3 15 1; 3 15 1])
%!error <FCN must be a valid function handle> ode45 (1, [0 25], [3 15 1])