Mercurial > octave
changeset 20631:00caf63edcdf
maint: Remove obsolete ODE options from odeset/odeget/ode45.
* odepkg_structure_check.m: Deleted file, replaced in code with
ode_struct_value_check()
* scripts/ode/module.mk: Remove odepkg_structure_check.m from build system.
* ode45.m: Just print_usage(), rather than displaying help, if incorrect
number of arguments supplied. Replace calls to odepkg_structure_check()
with calls to ode_struct_value_check(). Use error ID
"Octave:invalid-input-arg" rather than "OdePkg:InvalidArgument".
Use local variables TimeStepNumber and TimeStepSize since these are no
longer a default part of an ODE options structure. Remove options that
are no longer a part of ODE options structure.
* odeget.m: Remove unused options. Allow unknown options through, but
display a warning.
* odeset.m: Remove unused options. Allow unknown options through, but
dispaly a warning. Use special for loop over structure to simplify code.
Add BIST tests for custom user-defined options and for input validation.
* ode_struct_value_check.m: Rewrite docstring. Return input ode_struct
as output so that function can replace odepk_structure_check. Use
meaningful input variable name ode_struct rather than arg. Use special
for loop syntax over structure to simplify code. Remove checking for
obsolete, non-Matlab options. Use error ID "Octave:invalid-input-arg" rather
than "OdePkg:InvalidArgument". Preface all error and warning messages
with the name of the calling function.
* AbsRel_Norm.m, integrate_const.m, integrate_n_steps.m: Use error ID
"Octave:invalid-input-arg" rather than "OdePkg:InvalidArgument".
* odepkg_event_handle.m: Rephrase sentence in docstring.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 16 Oct 2015 10:31:02 -0700 |
parents | ab705b42cfd8 |
children | 80e630b37ba1 |
files | scripts/ode/module.mk scripts/ode/ode45.m scripts/ode/odeget.m scripts/ode/odeset.m scripts/ode/private/AbsRel_Norm.m scripts/ode/private/integrate_const.m scripts/ode/private/integrate_n_steps.m scripts/ode/private/ode_struct_value_check.m scripts/ode/private/odepkg_event_handle.m scripts/ode/private/odepkg_structure_check.m |
diffstat | 10 files changed, 216 insertions(+), 860 deletions(-) [+] |
line wrap: on
line diff
--- a/scripts/ode/module.mk Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/module.mk Fri Oct 16 10:31:02 2015 -0700 @@ -9,7 +9,6 @@ scripts/ode/private/integrate_n_steps.m \ scripts/ode/private/kahan.m \ scripts/ode/private/odepkg_event_handle.m \ - scripts/ode/private/odepkg_structure_check.m \ scripts/ode/private/ode_rk_interpolate.m \ scripts/ode/private/ode_struct_value_check.m \ scripts/ode/private/runge_kutta_45_dorpri.m \
--- a/scripts/ode/ode45.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/ode45.m Fri Oct 16 10:31:02 2015 -0700 @@ -23,7 +23,7 @@ ## @deftypefnx {Function File} {[@var{t}, @var{y}] =} ode45 (@var{fun}, @var{trange}, @var{init}, @var{opt}) ## @deftypefnx {Function File} {[@var{t}, @var{y}] =} ode45 (@dots{}, @var{par1}, @var{par2}, @dots{}) ## @deftypefnx {Function File} {[@var{t}, @var{y}, @var{xe}, @var{ye}, @var{ie}] =} ode45 (@dots{}) -## @deftypefnx {Function File} {@var{sol} =} ode45 (@var{fun}, @var{trange}, @var{init}, @dots{}) +## @deftypefnx {Function File} {@var{solution} =} ode45 (@var{fun}, @var{trange}, @var{init}, @dots{}) ## ## Solve a set of non-stiff Ordinary Differential Equations (non-stiff ODEs) ## with the well known explicit Dormand-Prince method of order 4. @@ -74,51 +74,44 @@ function varargout = ode45 (vfun, vtrange, vinit, varargin) - vorder = 5; # runge_kutta_45_dorpri uses local extrapolation - vsolver = "ode45"; - - ## FIXME: Octave core function's usually do print_usage () - ## rather than displaying full help string when improperly called. - if (nargin == 0) # Check number and types of all input arguments - help (vsolver); - error ("OdePkg:InvalidArgument", - "number of input arguments must be greater than zero"); - endif - if (nargin < 3) print_usage (); endif + vorder = 5; # runge_kutta_45_dorpri uses local extrapolation + vsolver = "ode45"; + if (nargin >= 4) if (! isstruct (varargin{1})) ## varargin{1:len} are parameters for vfun - vodeoptions = odeset; + vodeoptions = odeset (); vodeoptions.vfunarguments = varargin; elseif (length (varargin) > 1) - ## varargin{1} is an OdePkg options structure vopt - vodeoptions = odepkg_structure_check (varargin{1}, "ode45"); + ## varargin{1} is an ODE options structure vopt + vodeoptions = ode_struct_value_check ("ode45", varargin{1}, "ode45"); vodeoptions.vfunarguments = {varargin{2:length(varargin)}}; else # if (isstruct (varargin{1})) - vodeoptions = odepkg_structure_check (varargin{1}, "ode45"); + vodeoptions = ode_struct_value_check ("ode45", varargin{1}, "ode45"); vodeoptions.vfunarguments = {}; endif else # nargin == 3 - vodeoptions = odeset; + vodeoptions = odeset (); vodeoptions.vfunarguments = {}; endif if (! isvector (vtrange) || ! isnumeric (vtrange)) - error ("OdePkg:InvalidArgument", + error ("Octave:invalid-input-arg", "second input argument must be a valid vector"); endif + TimeStepNumber = odeget (vodeoptions, "TimeStepNumber", [], "fast"); + TimeStepSize = odeget (vodeoptions, "TimeStepSize", [], "fast"); if (length (vtrange) < 2 - && (isempty (vodeoptions.TimeStepSize) - || isempty (vodeoptions.TimeStepNumber))) - error ("OdePkg:InvalidArgument", + && (isempty (TimeStepSize) || isempty (TimeStepNumber))) + error ("Octave:invalid-input-arg", "second input argument must be a valid vector"); elseif (vtrange(2) == vtrange(1)) - error ("OdePkg:InvalidArgument", + error ("Octave:invalid-input-arg", "second input argument must be a valid vector"); else vodeoptions.vdirection = sign (vtrange(2) - vtrange(1)); @@ -126,63 +119,66 @@ vtrange = vtrange(:); if (! isvector (vinit) || ! isnumeric (vinit)) - error ("OdePkg:InvalidArgument", + error ("Octave:invalid-input-arg", "third input argument must be a valid numerical value"); endif vinit = vinit(:); if (ischar (vfun)) - try; vfun = str2func (vfun); catch; warning (lasterr); end_try_catch + try + vfun = str2func (vfun); + catch + warning (lasterr); + end_try_catch endif - if (! (isa (vfun, "function_handle"))) - error ("OdePkg:InvalidArgument", + if (! isa (vfun, "function_handle")) + error ("Octave:invalid-input-arg", "first input argument must be a valid function handle"); endif ## Start preprocessing, have a look which options are set in vodeoptions, ## check if an invalid or unused option is set - if (isempty (vodeoptions.TimeStepNumber) - && isempty (vodeoptions.TimeStepSize)) + if (isempty (TimeStepNumber) && isempty (TimeStepSize)) integrate_func = "adaptive"; vodeoptions.vstepsizefixed = false; - elseif (! isempty (vodeoptions.TimeStepNumber) - && ! isempty (vodeoptions.TimeStepSize)) + elseif (! isempty (TimeStepNumber) && ! isempty (TimeStepSize)) integrate_func = "n_steps"; vodeoptions.vstepsizefixed = true; - if (sign (vodeoptions.TimeStepSize) != vodeoptions.vdirection) - warning ("OdePkg:InvalidArgument", + if (sign (TimeStepSize) != vodeoptions.vdirection) + warning ("Octave:invalid-input-arg", "option 'TimeStepSize' has a wrong sign", "it will be corrected automatically"); - vodeoptions.TimeStepSize = (-1)*vodeoptions.TimeStepSize; + TimeStepSize = -TimeStepSize; endif - elseif (isempty (vodeoptions.TimeStepNumber) - && ! isempty (vodeoptions.TimeStepSize)) + elseif (isempty (TimeStepNumber) && ! isempty (TimeStepSize)) integrate_func = "const"; vodeoptions.vstepsizefixed = true; - if (sign (vodeoptions.TimeStepSize) != vodeoptions.vdirection) - warning ("OdePkg:InvalidArgument", + if (sign (TimeStepSize) != vodeoptions.vdirection) + warning ("Octave:invalid-input-arg", "option 'TimeStepSize' has a wrong sign", "it will be corrected automatically"); - vodeoptions.TimeStepSize = (-1)*vodeoptions.TimeStepSize; + TimeStepSize = -TimeStepSize; endif else - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "assuming an adaptive integrate function"); integrate_func = "adaptive"; endif + vodeoptions.TimeStepSize = TimeStepSize; + vodeoptions.TimeStepNumber = TimeStepNumber; ## Get the default options that can be set with "odeset" temporarily - vodetemp = odeset; + vodetemp = odeset (); ## Implementation of the option RelTol has been finished. ## This option can be set by the user to another value than default value. if (isempty (vodeoptions.RelTol) && ! vodeoptions.vstepsizefixed) vodeoptions.RelTol = 1e-3; - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'RelTol' not set, new value %f is used", vodeoptions.RelTol); elseif (! isempty (vodeoptions.RelTol) && vodeoptions.vstepsizefixed) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'RelTol' will be ignored if fixed time stamps are given"); endif @@ -190,11 +186,11 @@ ## This option can be set by the user to another value than default value. if (isempty (vodeoptions.AbsTol) && ! vodeoptions.vstepsizefixed) vodeoptions.AbsTol = 1e-6; - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'AbsTol' not set, new value %f is used", vodeoptions.AbsTol); elseif (! isempty (vodeoptions.AbsTol) && vodeoptions.vstepsizefixed) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'AbsTol' will be ignored if fixed time stamps are given"); else vodeoptions.AbsTol = vodeoptions.AbsTol(:); # Create column vector @@ -215,7 +211,7 @@ vodeoptions.vhavenonnegative = true; else vodeoptions.vhavenonnegative = false; - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'NonNegative' will be ignored if mass matrix is set"); endif else @@ -241,12 +237,6 @@ vodeoptions.vhaveoutputselection = false; endif - ## "OutputSave" option will be ignored. - if (! isempty (vodeoptions.OutputSave)) - warning ("OdePkg:InvalidArgument", - "Option 'OutputSave' will be ignored."); - endif - ## Implementation of the option Refine has been finished. ## This option can be set by the user to another value than default value. if (vodeoptions.Refine > 0) @@ -267,18 +257,18 @@ vodeoptions.AbsTol, vodeoptions.RelTol, vodeoptions.vnormcontrol); - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'InitialStep' not set, estimated value %f is used", vodeoptions.InitialStep); - elseif(isempty (vodeoptions.InitialStep)) - vodeoptions.InitialStep = odeget (vodeoptions, "TimeStepSize"); + elseif (isempty (vodeoptions.InitialStep)) + vodeoptions.InitialStep = TimeStepSize; endif ## Implementation of the option MaxStep has been finished. This option ## can be set by the user to another value than default value. if (isempty (vodeoptions.MaxStep) && ! vodeoptions.vstepsizefixed) vodeoptions.MaxStep = abs (vtrange(end) - vtrange(1)) / 10; - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "Option 'MaxStep' not set, new value %f is used", vodeoptions.MaxStep); endif @@ -295,30 +285,20 @@ ## by this solver because this solver uses an explicit Runge-Kutta method ## and therefore no Jacobian calculation is necessary. if (! isequal (vodeoptions.Jacobian, vodetemp.Jacobian)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'Jacobian' will be ignored by this solver"); endif if (! isequal (vodeoptions.JPattern, vodetemp.JPattern)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'JPattern' will be ignored by this solver"); endif if (! isequal (vodeoptions.Vectorized, vodetemp.Vectorized)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'Vectorized' will be ignored by this solver"); endif - if (! isequal (vodeoptions.NewtonTol, vodetemp.NewtonTol)) - warning ("OdePkg:InvalidArgument", - "option 'NewtonTol' will be ignored by this solver"); - endif - - if (! isequal (vodeoptions.MaxNewtonIterations,vodetemp.MaxNewtonIterations)) - warning ("OdePkg:InvalidArgument", - "option 'MaxNewtonIterations' will be ignored by this solver"); - endif - ## Implementation of the option Mass has been finished. ## This option can be set by the user to another value than default value. if (! isempty (vodeoptions.Mass) && isnumeric (vodeoptions.Mass)) @@ -341,27 +321,27 @@ ## Other options that are not used by this solver. ## Print a warning message to tell the user that the option is ignored. if (! isequal (vodeoptions.MvPattern, vodetemp.MvPattern)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'MvPattern' will be ignored by this solver"); endif if (! isequal (vodeoptions.MassSingular, vodetemp.MassSingular)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'MassSingular' will be ignored by this solver"); endif if (! isequal (vodeoptions.InitialSlope, vodetemp.InitialSlope)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'InitialSlope' will be ignored by this solver"); endif if (! isequal (vodeoptions.MaxOrder, vodetemp.MaxOrder)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'MaxOrder' will be ignored by this solver"); endif if (! isequal (vodeoptions.BDF, vodetemp.BDF)) - warning ("OdePkg:InvalidArgument", + warning ("Octave:invalid-input-arg", "option 'BDF' will be ignored by this solver"); endif @@ -387,12 +367,11 @@ case "n_steps" solution = integrate_n_steps (@runge_kutta_45_dorpri, vfun, vtrange(1), vinit, - vodeoptions.TimeStepSize, - vodeoptions.TimeStepNumber, SubOpts); + TimeStepSize, TimeStepNumber, SubOpts); case "const" solution = integrate_const (@runge_kutta_45_dorpri, vfun, vtrange, vinit, - vodeoptions.TimeStepSize, SubOpts); + TimeStepSize, SubOpts); endswitch ## Postprocessing, do whatever when terminating integration algorithm @@ -508,7 +487,7 @@ %! ## Turn off output of warning messages for all tests, turn them on %! ## again if the last test is called %!error # ouput argument -%! warning ("off", "OdePkg:InvalidArgument"); +%! warning ("off", "Octave:invalid-input-arg"); %! B = ode45 (1, [0 25], [3 15 1]); %!error # input argument number one %! [vt, vy] = ode45 (1, [0 25], [3 15 1]); @@ -532,11 +511,11 @@ %!test # extra input arguments passed through %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], 12, 13, "KL"); %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); -%!test # empty OdePkg structure *but* extra input arguments +%!test # empty ODEOPT structure *but* extra input arguments %! vopt = odeset; %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt, 12, 13, "KL"); %! assert ([vt(end), vy(end,:)], [2, fref], 1e-2); -%!error # strange OdePkg structure +%!error # strange ODEOPT structure %! vopt = struct ("foo", 1); %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt); %!test # Solve vdp in fixed step sizes @@ -639,8 +618,6 @@ %! %!## test for JPattern option is missing %!## test for Vectorized option is missing -%!## test for NewtonTol option is missing -%!## test for MaxNewtonIterations option is missing %! %!test # Mass option as function %! vopt = odeset ("Mass", @fmas); @@ -666,9 +643,10 @@ %! vopt = odeset ("BDF", "on"); %! [vt, vy] = ode45 (@fpol, [0 2], [2 0], vopt); %! assert ([vt(end), vy(end,:)], [2, fref], 1e-3); -%!test # +%!test %!## test for MvPattern option is missing %!## test for InitialSlope option is missing %!## test for MaxOrder option is missing %! -%! warning ("on", "OdePkg:InvalidArgument"); +%! warning ("on", "Octave:invalid-input-arg"); +
--- a/scripts/ode/odeget.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/odeget.m Fri Oct 16 10:31:02 2015 -0700 @@ -26,7 +26,7 @@ ## ## If called with two input arguments and the first input argument @var{ode_opt} ## is an ODE option structure and the second input argument @var{field} is a -## string specifying an option name then return the option value @var{val} +## string specifying an option name, then return the option value @var{val} ## corresponding to to @var{field} from @var{ode_opt}. ## ## If called called with an optional third input argument, and @var{field} is @@ -43,7 +43,7 @@ print_usage (); endif - ## Shortcut for empty options structures + ## Shortcut for empty option structures if (isempty (ode_opt)) if (nargin < 3) val = []; @@ -80,21 +80,16 @@ return; endif - ## Check if the given struct is a valid OdePkg struct - ode_struct_value_check (ode_opt); + ## Check if the given struct is a valid ODEOPT struct + ode_struct_value_check ("odeget", ode_opt); - ## Define all the possible OdePkg fields - persistent options = {"AbsTol"; "Algorithm"; "BDF"; "Choice"; "Eta"; "Events"; - "Explicit"; "InexactSolver"; "InitialSlope"; + ## Define all the possible ODEOPT fields + persistent options = {"AbsTol"; "BDF"; "Events"; "InitialSlope"; "InitialStep"; "Jacobian"; "JConstant"; "JPattern"; - "Mass"; "MassConstant"; "MassSingular"; - "MaxNewtonIterations"; "MaxOrder"; "MaxStep"; - "MStateDependence"; "MvPattern"; "NewtonTol"; - "NonNegative"; "NormControl"; "OutputFcn"; "OutputSave"; - "OutputSel"; "PolynomialDegree"; "QuadratureOrder"; - "Refine"; "RelTol"; "Restart"; "Stats"; - "TimeStepNumber"; "TimeStepSize"; "UseJacobian"; - "Vectorized"}; + "Mass"; "MassConstant"; "MassSingular"; "MaxOrder"; + "MaxStep"; "MStateDependence"; "MvPattern"; + "NonNegative"; "NormControl"; "OutputFcn"; "OutputSel"; + "Refine"; "RelTol"; "Stats"; "Vectorized"}; exactmatch = true; match = find (strcmpi (field, options)); @@ -104,15 +99,14 @@ endif if (isempty (match)) - if (nargin == 2) - error ("odeget: invalid property '%s'", field); - else - ## FIXME: Should we warn, but complete the action, or just error out? - warning ("odeget:InvalidArgument", - "odeget: invalid property '%s'. Using supplied default value.", - field); + ## Possibly a custom user-defined option + try + val = ode_opt.(field); + catch + warning ("Octave:invalid-input-arg", + "odeget: no field '%s' exists in ODE_OPT\n", field); val = default; - endif + end_try_catch elseif (numel (match) == 1) if (! exactmatch) warning ("odeget:NoExactMatching", @@ -135,7 +129,7 @@ %!demo -%! # Return the manually changed value RelTol of the OdePkg options +%! # Return the manually changed value RelTol of the ODE options %! # structure A. If RelTol wouldn't have been changed then an %! # empty matrix value would have been returned. %! @@ -157,8 +151,7 @@ %!error odeget (1,2,3,4,5) %!error <ODE_OPT must be a valid ODE_STRUCT> odeget (1, "opt1") %!error <FIELD must be a string> odeget (struct ("opt1", 1), 1) -%!error <invalid property 'foo'> odeget (struct ("opt1", 1), "foo") -%!warning <Using supplied default value> odeget (struct ("opt1", 1), "foo", 3); +%!warning <no field 'foo' exists> odeget (struct ("opt1", 1), "foo"); %!warning <no exact match for 'Rel'. Assuming 'RelTol'> odeget (struct ("RelTol", 1), "Rel"); %!error <Possible fields found: InitialSlope, InitialStep> odeget (odeset (), "Initial")
--- a/scripts/ode/odeset.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/odeset.m Fri Oct 16 10:31:02 2015 -0700 @@ -53,18 +53,13 @@ return; endif - ## Column vector of all possible OdePkg options - persistent options = {"AbsTol"; "Algorithm"; "BDF"; "Choice"; "Eta"; "Events"; - "Explicit"; "InexactSolver"; "InitialSlope"; + ## Column vector of all possible ODE options + persistent options = {"AbsTol"; "BDF"; "Events"; "InitialSlope"; "InitialStep"; "Jacobian"; "JConstant"; "JPattern"; - "Mass"; "MassConstant"; "MassSingular"; - "MaxNewtonIterations"; "MaxOrder"; "MaxStep"; - "MStateDependence"; "MvPattern"; "NewtonTol"; - "NonNegative"; "NormControl"; "OutputFcn"; "OutputSave"; - "OutputSel"; "PolynomialDegree"; "QuadratureOrder"; - "Refine"; "RelTol"; "Restart"; "Stats"; - "TimeStepNumber"; "TimeStepSize"; "UseJacobian"; - "Vectorized"}; + "Mass"; "MassConstant"; "MassSingular"; "MaxOrder"; + "MaxStep"; "MStateDependence"; "MvPattern"; + "NonNegative"; "NormControl"; "OutputFcn"; "OutputSel"; + "Refine"; "RelTol"; "Stats"; "Vectorized"}; ## initialize output odestruct = cell2struct (cell (numel (options), 1), options); @@ -75,13 +70,9 @@ if (isstruct (varargin{1})) oldstruct = varargin{1}; - ode_struct_value_check (oldstruct); - - oldstruct_fldnames = (fieldnames (oldstruct)).'; ## Copy oldstruct values into output odestruct - for fldname = oldstruct_fldnames - name = lower (fldname{1}); + for [val, name] = oldstruct exactmatch = true; match = find (strcmpi (name, options)); @@ -91,18 +82,25 @@ endif if (isempty (match)) - error ("odeset: invalid property '%s'", fldname{1}); + odestruct.(name) = val; elseif (numel (match) == 1) if (! exactmatch) warning ("odeset:NoExactMatching", "no exact match for '%s'. Assuming '%s'.", name, options{match}); endif - odestruct.(options{match}) = oldstruct.(fldname{1}); + odestruct.(options{match}) = val; else error ("odeset: no exact match for '%s'. Possible fields found: %s.", name, strjoin (options(match), ", ")); endif + + if (nargin == 1) + ## Check if all changes have resulted in a valid ODEOPT struct + ode_struct_value_check ("odeset", odestruct); + return; + endif + endfor ## At this point, odestruct has been initialized with default values, @@ -110,13 +108,9 @@ if (nargin == 2 && isstruct (varargin{2})) newstruct = varargin{2}; - ode_struct_value_check (newstruct); - - newstruct_fldnames = (fieldnames (newstruct)).'; ## Update the first struct with the values from the second one - for fldname = newstruct_fldnames - name = lower (fldname{1}); + for [val, name] = newstruct exactmatch = true; match = find (strcmpi (name, options)); @@ -126,21 +120,22 @@ endif if (isempty (match)) - error ("odeset: invalid property '%s'", fldname{1}); + odestruct.(name) = val; elseif (numel (match) == 1) if (! exactmatch) warning ("odeset:NoExactMatching", "no exact match for '%s'. Assuming '%s'.", name, options{match}); endif - odestruct.(options{match}) = newstruct.(fldname{1}); + odestruct.(options{match}) = val; else error ("odeset: no exact match for '%s'. Possible fields found: %s.", name, strjoin (options(match), ", ")); endif endfor - ## Done copying newstruct to oldstruct + ## Check if all changes have resulted in a valid ODEOPT struct + ode_struct_value_check ("odeset", odestruct); return; endif @@ -154,7 +149,7 @@ ## Write new field/value pairs into odestruct for i = 2:2:nargin - name = lower (varargin{i}); + name = varargin{i}; exactmatch = true; match = find (strcmpi (name, options)); @@ -164,7 +159,7 @@ endif if (isempty (match)) - error ("odeset: invalid property '%s'", varargin{i}); + odestruct.(name) = varargin{i+1}; elseif (numel (match) == 1) if (! exactmatch) warning ("odeset:NoExactMatching", @@ -178,8 +173,8 @@ endif endfor - ## Check if all changes have resulted in a valid OdePkg struct - ode_struct_value_check (odestruct); + ## Check if all changes have resulted in a valid ODEOPT struct + ode_struct_value_check ("odeset", odestruct); else ## First input argument was not a struct, must be field/value pairs @@ -190,7 +185,7 @@ endif for i = 1:2:nargin - name = lower (varargin{i}); + name = varargin{i}; exactmatch = true; match = find (strcmpi (name, options)); @@ -200,7 +195,7 @@ endif if (isempty (match)) - error ("odeset: invalid property '%s'", varargin{i}); + odestruct.(name) = varargin{i+1}; elseif (numel (match) == 1) if (! exactmatch) warning ("odeset:NoExactMatching", @@ -214,27 +209,22 @@ endif endfor - ## Check if all changes have resulted in a valid OdePkg struct - ode_struct_value_check (odestruct); + ## Check if all changes have resulted in a valid ODEOPT struct + ode_struct_value_check ("odeset", odestruct); endif endfunction -## function useful to print all the possible options +## function to print all possible options function print_options () disp ("List of all possible ODE solver options."); disp ("Default values are in square brackets."); disp (""); disp (" AbsTol: scalar or vector, >0, [1e-6]"); - disp (" Algorithm: string, {['gmres'], 'pcg', 'bicgstab'}"); disp (" BDF: binary, {'on', ['off']}"); - disp (" Choice: switch, {[1], 2}"); - disp (" Eta: scalar, >=0, <1, [0.5]"); disp (" Events: function_handle, []"); - disp (" Explicit: binary, {'yes', ['no']}"); - disp (" InexactSolver: string, {'inexact_newton', 'fsolve', []}"); disp (" InitialSlope: vector, []"); disp (" InitialStep: scalar, >0, []"); disp (" Jacobian: matrix or function_handle, []"); @@ -243,56 +233,46 @@ disp (" Mass: matrix or function_handle, []"); disp (" MassConstant: binary, {'on', ['off']}"); disp (" MassSingular: switch, {'yes', ['maybe'], 'no'}"); - disp ("MaxNewtonIterations: scalar, integer, >0, [1e3]"); disp (" MaxOrder: switch, {1, 2, 3, 4, [5]}"); disp (" MaxStep: scalar, >0, []"); disp (" MStateDependence: switch, {'none', ['weak'], 'strong'}"); disp (" MvPattern: sparse matrix, []"); - disp (" NewtonTol: scalar or vector, >0, []"); disp (" NonNegative: vector of integers, []"); disp (" NormControl: binary, {'on', ['off']}"); disp (" OutputFcn: function_handle, []"); - disp (" OutputSave: scalar, integer, >0, []"); disp (" OutputSel: scalar or vector, []"); - disp (" PolynomialDegree: scalar, integer, >0, []"); - disp (" QuadratureOrder: scalar, integer, >0, []"); disp (" Refine: scalar, integer, >0, []"); disp (" RelTol: scalar, >0, [1e-3]"); - disp (" Restart: scalar, integer, >0, [20]"); disp (" Stats: binary, {'on', ['off']}"); - disp (" TimeStepNumber: scalar, integer, >0, []"); - disp (" TimeStepSize: scalar, >0, []"); - disp (" UseJacobian: binary, {'yes', ['no']}"); disp (" Vectorized: binary, {'on', ['off']}"); endfunction %!demo -%! # A new OdePkg options structure with default values is created. +%! # A new ODE options structure with default values is created. %! %! odeoptA = odeset (); %!demo -%! # A new OdePkg options structure with manually set options +%! # 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 OdePkg options structure is created from odeoptB with +%! # 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"); -## All tests that are needed to check if a correct resp. valid option -## has been set are implemented in ode_struct_value_check.m. +## All tests that are needed to check if a valid option has been set are +## implemented in ode_struct_value_check.m %!test %! odeoptA = odeset (); %! assert (isstruct (odeoptA)); -%! fields = fieldnames (odeoptA); -%! assert (numel (fields), 37); +%! assert (numel (fieldnames (odeoptA)), 23); %! assert (all (structfun ("isempty", odeoptA))); %!shared odeoptB, odeoptC @@ -310,3 +290,21 @@ %! odeoptD = odeset (odeoptB, odeoptC); %! assert (odeoptD, odeoptC); +## Test custom user-defined option +%!test +%! wstate = warning ("off", "Octave:invalid-input-arg"); +%! unwind_protect +%! odeopt = odeset ("NewtonTol", 3); +%! assert (odeopt.NewtonTol, 3); +%! unwind_protect_cleanup +%! warning (wstate); +%! end_unwind_protect + +## Test input validation +%!error <FIELD/VALUE arguments must occur in pairs> odeset ("opt1") +%!error <FIELD names must be strings> odeset (1, 1) +%!error <FIELD/VALUE arguments must occur in pairs> odeset (odeset (), "opt1") +%!error <FIELD names must be strings> odeset (odeset (), 1, 1) +%!warning <no exact match for 'Rel'. Assuming 'RelTol'> odeset ("Rel", 1); +%!error <Possible fields found: InitialSlope, InitialStep> odeset ("Initial", 1) +
--- a/scripts/ode/private/AbsRel_Norm.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/private/AbsRel_Norm.m Fri Oct 16 10:31:02 2015 -0700 @@ -26,12 +26,12 @@ endif if (length (x_old) != n || length (y) != n) - error ("OdePkg:InvalidArgument", "invalid dimensions of input arguments"); + error ("Octave:invalid-input-arg", "invalid dimensions of input arguments"); endif if ((length (AbsTol) != 1 && length (AbsTol) != n) || (length (RelTol) != 1 && length (RelTol) != n)) - error ("OdePkg:InvalidArgument", "invalid dimensions of input arguments"); + error ("Octave:invalid-input-arg", "invalid dimensions of input arguments"); endif sc = AbsTol + max (abs (x), abs (x_old)) .* RelTol;
--- a/scripts/ode/private/integrate_const.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/private/integrate_const.m Fri Oct 16 10:31:02 2015 -0700 @@ -71,7 +71,7 @@ vdirection = odeget (options, "vdirection", [], "fast"); if (sign (dt) != vdirection) - error ("OdePkg:InvalidArgument", + error ("Octave:invalid-input-arg", "option 'InitialStep' has a wrong sign"); endif
--- a/scripts/ode/private/integrate_n_steps.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/private/integrate_n_steps.m Fri Oct 16 10:31:02 2015 -0700 @@ -72,7 +72,7 @@ vdirection = odeget (options, "vdirection", [], "fast"); if (sign (dt) != vdirection) - error ("OdePkg:InvalidArgument", "option 'InitialStep' has a wrong sign"); + error ("Octave:invalid-input-arg", "option 'InitialStep' has a wrong sign"); endif comp = 0.0;
--- a/scripts/ode/private/ode_struct_value_check.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/private/ode_struct_value_check.m Fri Oct 16 10:31:02 2015 -0700 @@ -18,31 +18,36 @@ ## <http://www.gnu.org/licenses/>. ## -*- texinfo -*- -## @deftypefn {Function File} {} ode_struct_value_check (@var{arg}) -## @deftypefnx {Function File} {} ode_struct_value_check (@var{arg}, @var{"solver"}) +## @deftypefn {Function File} {} ode_struct_value_check (@var{"caller"}, @var{ode_struct}) +## @deftypefnx {Function File} {} ode_struct_value_check (@var{"caller"), @var{ode_struct}, @var{"solver"}) +## @deftypefnx {Function File} {@var{ode_struct} =} ode_struct_value_check (@dots{}) +## +## Validate the fields and values in the ODE options structure @var{ode_struct}. +## +## The first argument @var{caller} is a string with the name of the calling +## function so that warning and error messages properly display the source +## of any problems. ## -## If this function is called with one input argument of type structure array -## then check the field names and the field values of the OdePkg structure -## @var{arg}. Optionally if this function is called with a second input -## argument @var{"solver"} of type string that specifies the name of a valid -## OdePkg solver then a higher level error detection is performed. The function -## does not modify any of the field names or field values but terminates with -## an error if an invalid option or value is found. +## The second argument @var{ode_struct} is a structure with fields and values +## that configure the ODE solvers (@pxref{XREFodeset,,odeset). +## +## The optional third argument @var{"solver"} is a string with the name of a +## specific ODE solver. This extra information can enable more extensive value +## validation for certain options. ## -## This function is an OdePkg internal helper function; Therefore, it should -## never be necessary for a user to call this function directly. +## The function does not modify any of the field names or field values, but +## terminates with an error if an invalid value is found. +## +## Normally the function is called with no output. However, the input struct +## is passed unmodified to the output for certain solvers which expect to +## receive the validated ODE structure returned. ## @end deftypefn ## ## @seealso{odeset, odeget} -function ode_struct_value_check (arg, solver = []) +function ode_struct = ode_struct_value_check (caller, ode_struct, solver = "") - fields = (fieldnames (arg)).'; - fields_nb = length (fields); - - for fldname = fields # Cycle over all fields - opt = fldname{1}; - val = arg.(opt); + for [val, opt] = ode_struct # Cycle over all fields switch (opt) @@ -50,67 +55,24 @@ if (! isempty (val)) if (! isnumeric (val) || ! isreal (val) || ! isvector (val) || any (val <= 0)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "Algorithm" - if (! isempty (val)) - if (! ischar (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "BDF" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "Choice" - if (! isempty (val)) - if (! isnumeric (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - elseif (val != 1 && val != 2) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "Eta" - if (! isempty (val)) - if (! isreal (val) || val < 0 || val >= 1) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "Events" if (! isempty (val)) if (! isa (val, "function_handle")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "Explicit" - if (! isempty (val)) - if (! strcmp (val, "yes") && ! strcmp (val, "no")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "InexactSolver" - if (! isempty (val)) - if (! ischar (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -118,8 +80,8 @@ if (! isempty (val)) if (! ischar (val) && (! isnumeric (val) || (! isvector (val) && ! isreal (val)))) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -127,8 +89,8 @@ if (! isempty (val)) if (! isnumeric (val) || ! isreal (val) || ! isscalar (val) || val <= 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -136,8 +98,8 @@ if (! isempty (val)) if (! isnumeric (val)) if (! isa (val, "function_handle") && ! iscell (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif endif @@ -145,16 +107,16 @@ case "JConstant" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "JPattern" if (! isempty (val)) if (! isnumeric (val) && ! isvector (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -162,33 +124,24 @@ if (! isempty (val)) if ((! isnumeric (val) || ! ismatrix (val)) && ! isa (val, "function_handle")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "MassConstant" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "MassSingular" if (! isempty (val)) if (! any (strcmp (val, {"yes", "no", "maybe"}))) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "MaxNewtonIterations" - if (! isempty (val)) - if (! isnumeric (val) - || val != fix (val) || val <= 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -196,40 +149,32 @@ if (! isempty (val)) if (! isnumeric (val) || val != fix (val) || val <= 0 || val >= 8) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "MaxStep" if (! isempty (val)) if (! isnumeric (val) || ! isscalar (val) || val <= 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "MStateDependence" if (! isempty (val)) if (! any (strcmp (val, {"none", "weak", "strong"}))) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "MvPattern" if (! isempty (val)) if (! isnumeric (val) && ! isvector (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "NewtonTol" - if (! isempty (val)) - if (! isnumeric (val) || ! isreal (val) || any (val <= 0)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -237,59 +182,32 @@ if (! isempty (val)) if (! isnumeric (val) || ! isvector (val) || any (val <= 0) || any (val != fix (val))) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "NormControl" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "OutputFcn" if (! isempty (val)) if (! isa (val, "function_handle")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "OutputSave" - if (! isempty (val)) - if (! isscalar (val) && val != Inf) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - elseif ((val != fix (val) || val <= 0) && val != Inf) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "OutputSel" if (! isempty (val)) if (! isnumeric (val) || ! isvector (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "PolynomialDegree" - if (! isempty (val)) - if (! isnumeric (val) || ! isvector (val) || any (val <= 0)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "QuadratureOrder" - if (! isempty (val)) - if (! isnumeric (val) || ! isvector (val) || any (val <= 0)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif @@ -297,77 +215,45 @@ if (! isempty (val)) if (! isnumeric (val) || ! isscalar (val) || val != fix (val) || val < 0 || val > 5) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "RelTol" if (! isempty (val)) if (! isnumeric (val) || ! isreal (val) || any (val <= 0)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif if (any (strcmp (solver, {"ode23", "ode23d", "ode45", "ode45d", "ode54", "ode54d", "ode78", "ode78d"}))) if (! isscalar (val)) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif endif - case "Restart" - if (! isempty (val)) - if (! isnumeric (val) || val != fix (val) || val <= 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - case "Stats" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "TimeStepNumber" - if (! isempty (val)) - if (val != fix (val) || val <= 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "TimeStepSize" - if (! isempty (val)) - if (! isreal (val) || val == 0) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); - endif - endif - - case "UseJacobian" - if (! isempty (val)) - if (! strcmp (val, "yes") && ! strcmp (val, "no")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif case "Vectorized" if (! isempty (val)) if (! strcmp (val, "on") && ! strcmp (val, "off")) - error ("OdePkg:InvalidArgument", - "invalid value assigned to field %s", opt); + error ("Octave:invalid-input-arg", + [caller ": invalid value assigned to field '%s'"], opt); endif endif otherwise - warning ("OdePkg:InvalidArgument", - "invalid field '%s' in ODE options", opt); + warning ("Octave:invalid-input-arg", + [caller ": unknown field '%s' in ODE options\n"], opt); endswitch endfor @@ -375,15 +261,15 @@ %!demo -%! # Return the checked OdePkg options structure that is created by +%! # Return the checked ODE options structure that is created by %! # the command odeset. %! %! ode_struct_value_check (odeset); %!demo -%! # Create the OdePkg options structure A with odeset and check it -%! # with odepkg_structure_check. This actually is unnecessary -%! # because odeset automatically calls odepkg_structure_check before +%! # Create the ODE options structure A with odeset and check it +%! # with ode_struct_value_check. This actually is unnecessary +%! # because odeset automatically calls ode_struct_value_check before %! # returning. %! %! A = odeset ();
--- a/scripts/ode/private/odepkg_event_handle.m Wed Oct 14 23:56:42 2015 -0400 +++ b/scripts/ode/private/odepkg_event_handle.m Fri Oct 16 10:31:02 2015 -0700 @@ -47,7 +47,7 @@ ## any type are given then pass these parameters through ## @command{odepkg_event_handle} to the event function. ## -## This function is an OdePkg internal helper function therefore it should +## This function is an ODE internal helper function therefore it should ## never be necessary that this function is called directly by a user. There ## is only little error detection implemented in this function file to ## achieve the highest performance.
--- a/scripts/ode/private/odepkg_structure_check.m Wed Oct 14 23:56:42 2015 -0400 +++ /dev/null Thu Jan 01 00:00:00 1970 +0000 @@ -1,498 +0,0 @@ -## Copyright (C) 2013 Roberto Porcu' <roberto.porcu@polimi.it> -## Copyright (C) 2006-2012 Thomas Treichl <treichl@users.sourceforge.net> -## -## 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 -## <http://www.gnu.org/licenses/>. - -## -*- texinfo -*- -## @deftypefn {Function File} {@var{newstruct} =} odepkg_structure_check (@var{oldstruct}) -## @deftypefnx {Function File} {@var{newstruct} =} odepkg_structure_check (@var{oldstruct}, @var{"solver"}) -## -## If this function is called with one input argument of type structure array -## then check the field names and the field values of the OdePkg structure -## @var{oldstruct} and return the structure as @var{newstruct} if no error is -## found. -## -## Optionally if this function is called with a second input argument -## @var{"solver"} of type string taht specifies the name of a valid OdePkg -## solver then a higher level error detection is performed. The function -## does not modify any of the field names or field values but terminates with -## an error if an invalid option or value is found. -## -## This function is an OdePkg internal helper function therefore it should -## never be necessary that this function is called directly by a user. There -## is only little error detection implemented in this function file to -## achieve the highest performance. -## -## Run examples with the command -## -## @example -## demo odepkg_structure_check -## @end example -## @end deftypefn -## -## @seealso{odepkg} - -function vret = odepkg_structure_check (varargin) - - ## Check the number of input arguments - if (nargin == 0) - help ("odepkg_structure_check"); - error ("OdePkg:InvalidArgument", - "Number of input arguments must be greater than zero"); - elseif (nargin > 2) - print_usage (); - elseif (nargin == 1 && isstruct (varargin{1})) - vret = varargin{1}; - vsol = ""; - vfld = fieldnames (vret); - vlen = length (vfld); - elseif (nargin == 2 && isstruct (varargin{1}) && ischar (varargin{2})) - vret = varargin{1}; - vsol = varargin{2}; - vfld = fieldnames (vret); - vlen = length (vfld); - endif - - for vcntarg = 1:vlen # Run through the number of given structure field names - - switch (vfld{vcntarg}) - - case "RelTol" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) && ... - isreal (vret.(vfld{vcntarg})) && ... - all (vret.(vfld{vcntarg}) > 0))) # "all" is a MatLab need - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - if (any (strcmp (vsol, {"ode23", "ode45", "ode54", "ode78", - "ode23d", "ode45d", "ode54d", "ode78d"}))) - if (! isscalar (vret.(vfld{vcntarg})) - && ! isempty (vret.(vfld{vcntarg}))) - error ("OdePkg:InvalidParameter", - 'Value of option "RelTol" must be a scalar'); - endif - endif - - case "AbsTol" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) && ... - isreal (vret.(vfld{vcntarg})) && ... - all (vret.(vfld{vcntarg}) > 0))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "NormControl" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "on") || ... - strcmp (vret.(vfld{vcntarg}), "off"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "NonNegative" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) - && isvector (vret.(vfld{vcntarg})))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "OutputFcn" - if (isempty (vret.(vfld{vcntarg})) || ... - isa (vret.(vfld{vcntarg}), "function_handle")) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "OutputSel" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) - && isvector (vret.(vfld{vcntarg}))) || ... - isscalar (vret.(vfld{vcntarg}))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "OutputSave" - if (isempty (vret.(vfld{vcntarg})) || ... - (isscalar (vret.(vfld{vcntarg})) && ... - mod (vret.(vfld{vcntarg}), 1) == 0 && ... - vret.(vfld{vcntarg}) > 0) || ... - vret.(vfld{vcntarg}) == Inf) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Refine" - if (isempty (vret.(vfld{vcntarg})) || ... - (isscalar (vret.(vfld{vcntarg})) && ... - mod (vret.(vfld{vcntarg}), 1) == 0 && ... - vret.(vfld{vcntarg}) >= 0 && ... - vret.(vfld{vcntarg}) <= 5)) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Stats" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "on") || ... - strcmp (vret.(vfld{vcntarg}), "off"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "InitialStep" - if (isempty (vret.(vfld{vcntarg})) || ... - (isscalar (vret.(vfld{vcntarg})) && ... - isreal (vret.(vfld{vcntarg})))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MaxStep" - if (isempty (vret.(vfld{vcntarg})) || ... - (isscalar (vret.(vfld{vcntarg})) && ... - vret.(vfld{vcntarg}) > 0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Events" - if (isempty (vret.(vfld{vcntarg})) || ... - isa (vret.(vfld{vcntarg}), "function_handle")) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Jacobian" - if (isempty (vret.(vfld{vcntarg})) || ... - isnumeric (vret.(vfld{vcntarg})) || ... - isa (vret.(vfld{vcntarg}), "function_handle") || ... - iscell (vret.(vfld{vcntarg}))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "JPattern" - if (isempty (vret.(vfld{vcntarg})) || ... - isvector (vret.(vfld{vcntarg})) || ... - isnumeric (vret.(vfld{vcntarg}))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Vectorized" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "on") || ... - strcmp (vret.(vfld{vcntarg}), "off"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Mass" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) || ... - isa (vret.(vfld{vcntarg}), "function_handle"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MStateDependence" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "none") || ... - strcmp (vret.(vfld{vcntarg}), "weak") || ... - strcmp (vret.(vfld{vcntarg}), "strong"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MvPattern" - if (isempty (vret.(vfld{vcntarg})) || ... - (isvector (vret.(vfld{vcntarg})) || ... - isnumeric (vret.(vfld{vcntarg})))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MassSingular" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "yes") || ... - strcmp (vret.(vfld{vcntarg}), "no") || ... - strcmp (vret.(vfld{vcntarg}), "maybe"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "InitialSlope" - if (isempty (vret.(vfld{vcntarg})) || ... - isvector (vret.(vfld{vcntarg}))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MaxOrder" - if (isempty (vret.(vfld{vcntarg})) || ... - (mod (vret.(vfld{vcntarg}), 1) == 0 && ... - vret.(vfld{vcntarg}) > 0 && ... - vret.(vfld{vcntarg}) < 8)) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "BDF" - if (isempty (vret.(vfld{vcntarg})) || ... - (strcmp (vret.(vfld{vcntarg}), "on") || ... - strcmp (vret.(vfld{vcntarg}), "off"))) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "NewtonTol" - if (isempty (vret.(vfld{vcntarg})) || ... - (isnumeric (vret.(vfld{vcntarg})) && ... - isreal (vret.(vfld{vcntarg})) && ... - all (vret.(vfld{vcntarg}) > 0))) # "all" is a MatLab need - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MaxNewtonIterations" - if (isempty (vret.(vfld{vcntarg})) || ... - (mod (vret.(vfld{vcntarg}), 1) == 0 && ... - vret.(vfld{vcntarg}) > 0)) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - ## new fields added - case "Algorithm" - if ( isempty (vret.(vfld{vcntarg})) || ischar (vret.(vfld{vcntarg})) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Choice" - if ( isempty (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) && (vret.(vfld{vcntarg})==1) - || vret.(vfld{vcntarg})==2 ) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Eta" - if ( isempty (vret.(vfld{vcntarg})) - || ( isreal (vret.(vfld{vcntarg})) - && vret.(vfld{vcntarg})>=0 && vret.(vfld{vcntarg})<1) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Explicit" - if ( isempty (vret.(vfld{vcntarg})) - || (ischar (vret.(vfld{vcntarg})) - && (strcmp (vret.(vfld{vcntarg}),"yes") - || strcmp (vret.(vfld{vcntarg}),"no"))) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "InexactSolver" - if ( isempty (vret.(vfld{vcntarg})) || ischar (vret.(vfld{vcntarg})) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "InitialSlope" - if ( isempty (vret.(vfld{vcntarg})) - || ( ischar (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) - && (isvector (vret.(vfld{vcntarg})) - || isreal (vret.(vfld{vcntarg}))))) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "JConstant" - if ( isempty (vret.(vfld{vcntarg})) - || (ischar (vret.(vfld{vcntarg})) - && (strcmp (vret.(vfld{vcntarg}),"yes") - || strcmp (vret.(vfld{vcntarg}),"no"))) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "MassConstant" - if ( isempty (vret.(vfld{vcntarg})) - || (strcmp (vret.(vfld{vcntarg}),"on") - || strcmp (vret.(vfld{vcntarg}),"off")) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "PolynomialDegree" - if ( isempty (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) - && mod(vret.(vfld{vcntarg}),1)==0 && vret.(vfld{vcntarg})>0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "QuadratureOrder" - if ( isempty (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) - && mod(vret.(vfld{vcntarg}),1)==0 && vret.(vfld{vcntarg})>0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "Restart" - if ( isempty (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) - && mod(vret.(vfld{vcntarg}),1)==0 && vret.(vfld{vcntarg})>0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "TimeStepNumber" - if ( isempty (vret.(vfld{vcntarg})) - || (isnumeric (vret.(vfld{vcntarg})) - && mod(vret.(vfld{vcntarg}),1)==0 && vret.(vfld{vcntarg})>0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "TimeStepSize" - if ( isempty (vret.(vfld{vcntarg})) - || ( isreal (vret.(vfld{vcntarg})) && vret.(vfld{vcntarg})!=0) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - case "UseJacobian" - if ( isempty (vret.(vfld{vcntarg})) - || (ischar (vret.(vfld{vcntarg})) - && (strcmp (vret.(vfld{vcntarg}),"yes") - || strcmp (vret.(vfld{vcntarg}),"no"))) ) - else - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s" or no valid parameter value', - vfld{vcntarg}); - endif - - otherwise - error ("OdePkg:InvalidParameter", - 'Unknown parameter name "%s"', - vfld{vcntarg}); - - endswitch - - endfor - -endfunction - - -%!demo -%! # Return the checked OdePkg options structure that is created by -%! # the command odeset. -%! -%! odepkg_structure_check (odeset); - -%!demo -%! # Create the OdePkg options structure A with odeset and check it -%! # with odepkg_structure_check. This actually is unnecessary -%! # because odeset automatically calls odepkg_structure_check before -%! # returning. -%! -%! A = odeset (); odepkg_structure_check (A); -