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
|
########################################################################
##
## 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{odestruct} =} odeset ()
## @deftypefnx {} {@var{odestruct} =} odeset (@var{"field1"}, @var{value1}, @var{"field2"}, @var{value2}, @dots{})
## @deftypefnx {} {@var{odestruct} =} odeset (@var{oldstruct}, @var{"field1"}, @var{value1}, @var{"field2"}, @var{value2}, @dots{})
## @deftypefnx {} {@var{odestruct} =} odeset (@var{oldstruct}, @var{newstruct})
## @deftypefnx {} {} odeset ()
##
## Create or modify an ODE options structure.
##
## When called with no input argument and one output argument, return a new ODE
## options structure that contains all possible fields initialized to their
## default values. If no output argument is requested, display a list of
## the common ODE solver options along with their default value.
##
## If called with name-value input argument pairs @var{"field1"},
## @var{"value1"}, @var{"field2"}, @var{"value2"}, @dots{} return a new
## ODE options structure with all the most common option fields
## initialized, @strong{and} set the values of the fields @var{"field1"},
## @var{"field2"}, @dots{} to the values @var{value1}, @var{value2},
## @dots{}.
##
## If called with an input structure @var{oldstruct} then overwrite the
## values of the options @var{"field1"}, @var{"field2"}, @dots{} with
## new values @var{value1}, @var{value2}, @dots{} and return the
## modified structure.
##
## When called with two input ODE options structures @var{oldstruct} and
## @var{newstruct} overwrite all values from the structure
## @var{oldstruct} with new values from the structure @var{newstruct}.
## Empty values in @var{newstruct} will not overwrite values in
## @var{oldstruct}.
##
## The most commonly used ODE options, which are always assigned a value
## by @code{odeset}, are the following:
##
## @table @asis
## @item @code{AbsTol}: positive scalar | vector, def. @code{1e-6}
## Absolute error tolerance.
##
## @item @code{BDF}: @{@qcode{"off"}@} | @qcode{"on"}
## Use BDF formulas in implicit multistep methods.
## @emph{Note}: This option is not yet implemented.
##
## @item @code{Events}: function_handle
## Event function. An event function must have the form
## @code{[value, isterminal, direction] = my_events_f (t, y)}
##
## @item @code{InitialSlope}: vector
## Consistent initial slope vector for DAE solvers.
##
## @item @code{InitialStep}: positive scalar
## Initial time step size.
##
## @item @code{Jacobian}: matrix | function_handle
## Jacobian matrix, specified as a constant matrix or a function of
## time and state.
##
## @item @code{JConstant}: @{@qcode{"off"}@} | @qcode{"on"}
## Specify whether the Jacobian is a constant matrix or depends on the
## state.
##
## @item @code{JPattern}: sparse matrix
## If the Jacobian matrix is sparse and non-constant but maintains a
## constant sparsity pattern, specify the sparsity pattern.
##
## @item @code{Mass}: matrix | function_handle
## Mass matrix, specified as a constant matrix or a function of
## time and state.
##
## @item @code{MassSingular}: @{@qcode{"maybe"}@} | @qcode{"yes"} | @qcode{"on"}
## Specify whether the mass matrix is singular.
##
## @item @code{MaxOrder}: @{@qcode{5}@} | @qcode{4} | @qcode{3} | @qcode{2} | @qcode{1}
## Maximum order of formula.
##
## @item @code{MaxStep}: positive scalar
## Maximum time step value.
##
## @item @code{MStateDependence}: @{@qcode{"weak"}@} | @qcode{"none"} | @qcode{"strong"}
## Specify whether the mass matrix depends on the state or only on time.
##
## @item @code{MvPattern}: sparse matrix
## If the mass matrix is sparse and non-constant but maintains a
## constant sparsity pattern, specify the sparsity pattern.
## @emph{Note}: This option is not yet implemented.
##
## @item @code{NonNegative}: scalar | vector
## Specify elements of the state vector that are expected to remain
## non-negative during the simulation.
##
## @item @code{NormControl}: @{@qcode{"off"}@} | @qcode{"on"}
## Control error relative to the 2-norm of the solution, rather than its
## absolute value.
##
## @item @code{OutputFcn}: function_handle
## Function to monitor the state during the simulation. For the form of
## the function to use @pxref{XREFodeplot,,@code{odeplot}}.
##
## @item @code{OutputSel}: scalar | vector
## Indices of elements of the state vector to be passed to the output
## monitoring function.
##
## @item @code{Refine}: positive scalar
## Specify whether output should be returned only at the end of each
## time step or also at intermediate time instances. The value should be
## a scalar indicating the number of equally spaced time points to use
## within each timestep at which to return output.
##
## @item @code{RelTol}: positive scalar
## Relative error tolerance.
##
## @item @code{Stats}: @{@qcode{"off"}@} | @qcode{"on"}
## Print solver statistics after simulation.
##
## @item @code{Vectorized}: @{@qcode{"off"}@} | @qcode{"on"}
## Specify whether @code{odefcn} can be passed multiple values of the
## state at once.
##
## @end table
##
## Field names that are not in the above list are also accepted and
## added to the result structure.
##
## @seealso{odeget}
## @end deftypefn
function odestruct = odeset (varargin)
persistent p;
if (isempty (p))
p = inputParser ();
p.addParameter ("AbsTol", []);
p.addParameter ("BDF", []);
p.addParameter ("Events", []);
p.addParameter ("InitialSlope", []);
p.addParameter ("InitialStep", []);
p.addParameter ("Jacobian", []);
p.addParameter ("JConstant", []);
p.addParameter ("JPattern", []);
p.addParameter ("Mass", []);
p.addParameter ("MassSingular", []);
p.addParameter ("MaxOrder", []);
p.addParameter ("MaxStep", []);
p.addParameter ("MStateDependence", []);
p.addParameter ("MvPattern", []);
p.addParameter ("NonNegative", []);
p.addParameter ("NormControl", []);
p.addParameter ("OutputFcn", []);
p.addParameter ("OutputSel", []);
p.addParameter ("Refine", []);
p.addParameter ("RelTol", []);
p.addParameter ("Stats", []);
p.addParameter ("Vectorized", []);
p.KeepUnmatched = true;
endif
if (nargin == 0 && nargout == 0)
print_options ();
else
p.parse (varargin{:});
odestruct = p.Results;
odestruct_extra = p.Unmatched;
xtra_fields = fieldnames (odestruct_extra);
if (! isempty (xtra_fields))
## Merge extra fields into existing odestruct
for fldname = sort (xtra_fields.')
fldname = fldname{1};
warning ("Octave:invalid-input-arg",
"odeset: unknown option \"%s\"\n", fldname);
odestruct.(fldname) = odestruct_extra.(fldname);
endfor
endif
endif
endfunction
## function to print all possible options
function print_options ()
disp ("List of the most common ODE solver options.");
disp ("Default values are in square brackets.");
disp ("");
disp (' AbsTol: scalar or vector, >0, [1e-6]');
disp (' BDF: binary, {["off"], "on"}');
disp (' Events: function_handle, []');
disp (' InitialSlope: vector, []');
disp (' InitialStep: scalar, >0, []');
disp (' Jacobian: matrix or function_handle, []');
disp (' JConstant: binary, {["off"], "on"}');
disp (' JPattern: sparse matrix, []');
disp (' Mass: matrix or function_handle, []');
disp (' MassSingular: switch, {["maybe"], "no", "yes"}');
disp (' MaxOrder: switch, {[5], 1, 2, 3, 4, }');
disp (' MaxStep: scalar, >0, []');
disp (' MStateDependence: switch, {["weak"], "none", "strong"}');
disp (' MvPattern: sparse matrix, []');
disp (' NonNegative: vector of integers, []');
disp (' NormControl: binary, {["off"], "on"}');
disp (' OutputFcn: function_handle, []');
disp (' OutputSel: scalar or vector, []');
disp (' Refine: scalar, integer, >0, []');
disp (' RelTol: scalar, >0, [1e-3]');
disp (' Stats: binary, {["off"], "on"}');
disp (' Vectorized: binary, {["off"], "on"}');
endfunction
%!demo
%! ## A new ODE options structure with default values is created.
%!
%! odeoptA = odeset ();
%!demo
%! ## A new ODE options structure with manually set options
%! ## for "AbsTol" and "RelTol" is created.
%!
%! odeoptB = odeset ("AbsTol", 1e-2, "RelTol", 1e-1);
%!demo
%! ## A new ODE options structure is created from odeoptB with
%! ## a modified value for option "NormControl".
%!
%! odeoptB = odeset ("AbsTol", 1e-2, "RelTol", 1e-1);
%! odeoptC = odeset (odeoptB, "NormControl", "on");
%!test
%! odeoptA = odeset ();
%! assert (isstruct (odeoptA));
%! assert (numfields (odeoptA), 22);
%! assert (all (structfun ("isempty", odeoptA)));
%!shared odeoptB, odeoptC
%!test
%! odeoptB = odeset ("ABSTOL", 1e-2, "reltol", 1e-1);
%! assert (odeoptB.AbsTol, 1e-2); # Check canonicalization of name
%! assert (odeoptB.RelTol, 1e-1);
%!test
%! odeoptC = odeset (odeoptB, "NormControl", "on");
%! assert (odeoptC.AbsTol, 1e-2); # check values from first struct copied
%! assert (odeoptC.NormControl, "on"); # check new values override old ones
%!test
%! odeoptD = odeset (odeoptB, odeoptC);
%! assert (odeoptD, odeoptC);
## Test custom user-defined option
%!test
%! warning ("off", "Octave:invalid-input-arg", "local");
%! odeopt = odeset ("NewtonTol", 3);
%! assert (odeopt.NewtonTol, 3);
## Test input validation
%!error <argument 'opt1' is not a declared parameter> odeset ("opt1")
%!error odeset (1, 1)
%!error <argument 'opt1' is not a declared parameter> odeset (odeset (), "opt1")
%!error odeset (odeset (), 1, 1)
%!error <'Re' matches more than one Parameter: 'Refine', 'RelTol'>
%! odeset ('Re', 1);
|