Mercurial > octave
changeset 21443:acd6e203031d
Alter BIST tests stop emitting warnings during runtests invocation.
* hsv2rgb.m, rgb2hsv.m, rgb2ntsc.m: Use intended variable of 'int8' rather
than 'uint16'.
* strread.m: Change %!test to %!warning to verify, but swallow, emitted
warning. Wrap test to < 80 characters. Place comments about bugs
above %!test block. Use %!assert blocks rather than %!test, followed
by assert. Add missing semicolons.
* ode23.m: Clean up file to follow Octave guidelines for indentation.
Use '#' for comments which follow code on a line.
Use 'warning ("off", ..., "local")' to temporarily disable warning messages
during BIST testing.
* ode45.m: Add %!error tests at end of file.
Use 'warning ("off", ..., "local")' to temporarily disable warning messages
during BIST testing.
* pcg.m: Change %!test to %!warning to verify, but swallow, emitted warning.
* hankel.m: Change %!test to %!warning to verify, but swallow, emitted warning.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 11 Mar 2016 17:12:17 -0800 |
parents | 78edcabed953 |
children | 74b07936fd05 |
files | scripts/image/hsv2rgb.m scripts/image/rgb2hsv.m scripts/image/rgb2ntsc.m scripts/io/strread.m scripts/ode/ode23.m scripts/ode/ode45.m scripts/sparse/pcg.m scripts/special-matrix/hankel.m |
diffstat | 8 files changed, 229 insertions(+), 219 deletions(-) [+] |
line wrap: on
line diff
--- a/scripts/image/hsv2rgb.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/image/hsv2rgb.m Fri Mar 11 17:12:17 2016 -0800 @@ -182,7 +182,7 @@ %! assert (size (rgb), [10 10 3]) %!test -%! rgb = hsv2rgb (randi ([-128 127], 10, 10, 3, "uint16")); +%! rgb = hsv2rgb (randi ([-128 127], 10, 10, 3, "int8")); %! assert (class (rgb), "double") %! assert (size (rgb), [10 10 3])
--- a/scripts/image/rgb2hsv.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/image/rgb2hsv.m Fri Mar 11 17:12:17 2016 -0800 @@ -158,7 +158,7 @@ %! assert (size (hsv), [10 10 3]) %!test -%! hsv = rgb2hsv (randi ([-128 127], 10, 10, 3, "uint16")); +%! hsv = rgb2hsv (randi ([-128 127], 10, 10, 3, "int8")); %! assert (class (hsv), "double") %! assert (size (hsv), [10 10 3]) @@ -172,3 +172,4 @@ %! assert (rgb2hsv (rgb_double), expected) %! assert (rgb2hsv (rgb_uint8), expected, 0.005) %! assert (rgb2hsv (single (rgb_double)), single (expected)) +
--- a/scripts/image/rgb2ntsc.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/image/rgb2ntsc.m Fri Mar 11 17:12:17 2016 -0800 @@ -136,7 +136,7 @@ %! assert (size (ntsc), [10 10 3]) %!test -%! ntsc = rgb2ntsc (randi ([-128 127], 10, 10, 3, "uint16")); +%! ntsc = rgb2ntsc (randi ([-128 127], 10, 10, 3, "int8")); %! assert (class (ntsc), "double") %! assert (size (ntsc), [10 10 3])
--- a/scripts/io/strread.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/io/strread.m Fri Mar 11 17:12:17 2016 -0800 @@ -870,23 +870,23 @@ %! assert (b, {"2"}); %!test -%! assert (strread ("Hello World! // this is comment", "%s",... -%! "commentstyle", "c++"), ... -%! {"Hello"; "World!"}); +%! assert (strread ("Hello World! // this is comment", "%s", +%! "commentstyle", "c++"), +%! {"Hello"; "World!"}); %! assert (strread ("Hello World! % this is comment", "%s",... -%! "commentstyle", "matlab"), ... -%! {"Hello"; "World!"}); +%! "commentstyle", "matlab"), ... +%! {"Hello"; "World!"}); %! assert (strread ("Hello World! # this is comment", "%s",... -%! "commentstyle", "shell"), ... -%! {"Hello"; "World!"}); +%! "commentstyle", "shell"), ... +%! {"Hello"; "World!"}); %!test %! str = sprintf ("Tom 100 miles/hr\nDick 90 miles/hr\nHarry 80 miles/hr"); %! fmt = "%s %f miles/hr"; %! c = cell (1, 2); %! [c{:}] = strread (str, fmt); -%! assert (c{1}, {"Tom"; "Dick"; "Harry"}) -%! assert (c{2}, [100; 90; 80]) +%! assert (c{1}, {"Tom"; "Dick"; "Harry"}); +%! assert (c{2}, [100; 90; 80]); %!test %! a = strread ("a b c, d e, , f", "%s", "delimiter", ","); @@ -904,15 +904,15 @@ %! assert (a, int32 (10)); %! assert (b, {"a"}); +## Bug #33536 %!test -%! ## Bug #33536 %! [a, b, c] = strread ("1,,2", "%s%s%s", "delimiter", ","); %! assert (a{1}, "1"); %! assert (b{1}, ""); %! assert (c{1}, "2"); +## Bug #33536 %!test -%! ## Bug #33536 %! a = strread ("[SomeText]", "[%s", "delimiter", "]"); %! assert (a{1}, "SomeText"); @@ -932,8 +932,8 @@ %! assert (c, [3; NaN]); %! assert (d, int32 ([0; 0])); +## Default format (= %f) %!test -%! ## Default format (= %f) %! [a, b, c] = strread ("0.12 0.234 0.3567"); %! assert (a, 0.12); %! assert (b, 0.234); @@ -943,14 +943,14 @@ %! [a, b] = strread ("0.41 8.24 3.57 6.24 9.27", "%f%f", 2, "delimiter", " "); %! assert (a, [0.41; 3.57]); +## TreatAsEmpty %!test -%! ## TreatAsEmpty %! [a, b, c, d] = strread ("1,2,3,NN,5,6\n", "%d%d%d%f", "delimiter", ",", "TreatAsEmpty", "NN"); %! assert (c, int32 ([3; 0])); %! assert (d, [NaN; NaN]); +## No delimiters at all besides EOL. Plain reading numbers & strings %!test -%! ## No delimiters at all besides EOL. Plain reading numbers & strings %! str = "Text1Text2Text\nText398Text4Text\nText57Text"; %! [a, b] = strread (str, "Text%dText%1sText"); %! assert (a, int32 ([1; 398; 57])); @@ -967,55 +967,52 @@ %! assert (d', [15, 25, 35]); ## Bug #44750 -%!test -%! assert (strread ('/home/foo/','%s','delimiter','/','MultipleDelimsAsOne',1), ... -%! {"home"; "foo"}); +%!assert (strread ('/home/foo/','%s','delimiter','/','MultipleDelimsAsOne',1), +%! {"home"; "foo"}) ## delimiter as sq_string and dq_string -%!test -%! assert (strread ("1\n2\n3", "%d", "delimiter", "\n"), -%! strread ("1\n2\n3", "%d", "delimiter", '\n')) +%!assert (strread ("1\n2\n3", "%d", "delimiter", "\n"), +%! strread ("1\n2\n3", "%d", "delimiter", '\n')) ## whitespace as sq_string and dq_string -%!test -%! assert (strread ("1\b2\r3\b4\t5", "%d", "whitespace", "\b\r\n\t"), -%! strread ("1\b2\r3\b4\t5", "%d", "whitespace", '\b\r\n\t')) +%!assert (strread ("1\b2\r3\b4\t5", "%d", "whitespace", "\b\r\n\t"), +%! strread ("1\b2\r3\b4\t5", "%d", "whitespace", '\b\r\n\t')) %!test %! str = "0.31 0.86 0.94\n 0.60 0.72 0.87"; %! fmt = "%f %f %f"; %! args = {"delimiter", " ", "endofline", "\n", "whitespace", " "}; %! [a, b, c] = strread (str, fmt, args{:}); -%! assert (a, [0.31; 0.60], 0.01) -%! assert (b, [0.86; 0.72], 0.01) -%! assert (c, [0.94; 0.87], 0.01) +%! assert (a, [0.31; 0.60], 0.01); +%! assert (b, [0.86; 0.72], 0.01); +%! assert (c, [0.94; 0.87], 0.01); %!test %! str = "0.31,0.86,0.94\n0.60,0.72,0.87"; %! fmt = "%f %f %f"; %! args = {"delimiter", ",", "endofline", "\n", "whitespace", " "}; %! [a, b, c] = strread (str, fmt, args{:}); -%! assert (a, [0.31; 0.60], 0.01) -%! assert (b, [0.86; 0.72], 0.01) -%! assert (c, [0.94; 0.87], 0.01) +%! assert (a, [0.31; 0.60], 0.01); +%! assert (b, [0.86; 0.72], 0.01); +%! assert (c, [0.94; 0.87], 0.01); %!test %! str = "0.31 0.86 0.94\n 0.60 0.72 0.87"; %! fmt = "%f %f %f"; %! args = {"delimiter", ",", "endofline", "\n", "whitespace", " "}; %! [a, b, c] = strread (str, fmt, args{:}); -%! assert (a, [0.31; 0.60], 0.01) -%! assert (b, [0.86; 0.72], 0.01) -%! assert (c, [0.94; 0.87], 0.01) +%! assert (a, [0.31; 0.60], 0.01); +%! assert (b, [0.86; 0.72], 0.01); +%! assert (c, [0.94; 0.87], 0.01); %!test %! str = "0.31, 0.86, 0.94\n 0.60, 0.72, 0.87"; %! fmt = "%f %f %f"; %! args = {"delimiter", ",", "endofline", "\n", "whitespace", " "}; %! [a, b, c] = strread (str, fmt, args{:}); -%! assert (a, [0.31; 0.60], 0.01) -%! assert (b, [0.86; 0.72], 0.01) -%! assert (c, [0.94; 0.87], 0.01) +%! assert (a, [0.31; 0.60], 0.01); +%! assert (b, [0.86; 0.72], 0.01); +%! assert (c, [0.94; 0.87], 0.01); %!test %! [a, b] = strread (["Empty 1" char(10)], "Empty%s %f"); @@ -1027,22 +1024,22 @@ %! assert (a, NaN); %! assert (b, NaN); +## Bug #35999 %!test -%! ## Bug #35999 %! [a, b, c] = strread ("", "%f"); %! assert (isempty (a)); %! assert (isempty (b)); %! assert (isempty (c)); +## Bug #37023 %!test -%! ## bug #37023 %! [a, b] = strread (" 1. 1 \n 2 3 \n", "%f %f", "endofline", "\n"); %! assert (a, [1; 2], 1e-15); %! assert (b, [1; 3], 1e-15); ## Test for no output arg (interactive use) -%!test -%! assert (strread (",2,,4\n5,,7,", "", "delimiter", ","), [NaN; 2; NaN; 4; 5; NaN; 7]); +%!assert (strread (",2,,4\n5,,7,", "", "delimiter", ","), +%! [NaN; 2; NaN; 4; 5; NaN; 7]) ## Test #1 bug #42609 %!test @@ -1073,7 +1070,6 @@ %! assert (c, [3; 6; 9]); ## Unsupported format specifiers -%!test %!error <format specifiers are not supported> strread ("a", "%c") %!error <format specifiers are not supported> strread ("a", "%*c %d") %!error <format specifiers are not supported> strread ("a", "%q") @@ -1090,12 +1086,10 @@ %!error <format specifiers are not supported> strread ("a", "%*u32 %d") ## Illegal format specifiers -%!test -%!error <no valid format conversion specifiers> strread ("1.0", "%z"); +%!error <no valid format conversion specifiers> strread ("1.0", "%z") ## Test for false positives in check for non-supported format specifiers -%!test -%! assert (strread ("Total: 32.5 % (of cm values)","Total: %f % (of cm values)"), 32.5, 1e-5); +%!assert (strread ("Total: 32.5 % (of cm values)","Total: %f % (of cm values)"), 32.5, 1e-5) ## Test various forms of string format specifiers (bug #45712) %!test @@ -1103,8 +1097,8 @@ %! [a, b, c, d] = strread (str, "%f %s %*s %3s %*3s %f", "delimiter", ":"); %! assert ({a, b, c, d}, {14, {"1 z"}, {"3 z"}, 11}); -## Allow cuddling %sliteral but warn it is ambiguous -%!test +## Allow cuddling %sliteral but warn that it is ambiguous +%!warning <Ambiguous '%s' specifier immediately before literal in column 1> %! [a, b] = strread ("abcxyz51\nxyz83\n##xyz101", "%s xyz %d"); %! assert (a([1 3]), {"abc"; "##"}); %! assert (isempty (a{2}), true);
--- a/scripts/ode/ode23.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/ode/ode23.m Fri Mar 11 17:12:17 2016 -0800 @@ -86,7 +86,7 @@ ## @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}] = ode23 (fvdp, [0 20], [2 0]); +## [@var{t},@var{y}] = ode23 (fvdp, [0, 20], [2, 0]); ## @end group ## @end example ## @seealso{odeset, odeget} @@ -127,7 +127,7 @@ odeopts.funarguments = {}; endif - if (! isvector (trange) || ! isnumeric (trange)) + if (! isnumeric (trange) || ! isvector (trange)) error ("Octave:invalid-input-arg", "ode23: TRANGE must be a numeric vector"); endif @@ -146,7 +146,7 @@ endif trange = trange(:); - if (! isvector (init) || ! isnumeric (init)) + if (! isnumeric (init) || ! isvector (init)) error ("Octave:invalid-input-arg", "ode23: INIT must be a numeric vector"); endif @@ -312,7 +312,7 @@ "ode23: option \"BDF\" is ignored by this solver\n"); endif - ## Starting the initialisation of the core solver ode23 + ## Starting the initialization of the core solver ode23 if (havemasshandle) # Handle only the dynamic mass matrix, if (massdependence) # constant mass matrices have already @@ -356,9 +356,9 @@ nsteps = solution.cntloop-2; # cntloop from 2..end nfailed = (solution.cntcycles-1)-nsteps+1; # cntcycl from 1..end nfevals = 3 * (solution.cntcycles-1); # number of ode evaluations - ndecomps = 0; # number of LU decompositions - npds = 0; # number of partial derivatives - nlinsols = 0; # no. of solutions of linear systems + ndecomps = 0; # number of LU decompositions + npds = 0; # number of partial derivatives + nlinsols = 0; # no. of solutions of linear systems ## Print cost statistics if no output argument is given if (nargout == 0) printf ("Number of successful steps: %d\n", nsteps); @@ -411,13 +411,13 @@ %!function ydot = fpol (t, y) # The Van der Pol %! ydot = [y(2); (1 - y(1)^2) * y(2) - y(1)]; %!endfunction -%!function ref = fref () # The computed reference sol +%!function ref = fref () # The computed reference sol %! ref = [0.32331666704577, -1.83297456798624]; %!endfunction -%!function jac = fjac (t, y, varargin) # its Jacobian +%!function jac = fjac (t, y, varargin) # its Jacobian %! jac = [0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2]; %!endfunction -%!function jac = fjcc (t, y, varargin) # sparse type +%!function jac = fjcc (t, y, varargin) # sparse type %! jac = sparse ([0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2]) %!endfunction %!function [val, trm, dir] = feve (t, y, varargin) @@ -449,140 +449,144 @@ %! endif %!endfunction %! -%! ## Turn off output of warning messages for all tests, -%! ## turn them on again after the last test is called -%!test -%! warning ("off", "Octave:invalid-input-arg", "local"); -%!error ## ouput argument -%! B = ode23 (1, [0 25], [3 15 1]); -%!error ## input argument number one -%! [t, y] = ode23 (1, [0 25], [3 15 1]); -%!error ## input argument number two -%! [t, y] = ode23 (@fpol, 1, [3 15 1]); -%!test ## two output arguments -%! [t, y] = ode23 (@fpol, [0 2], [2 0]); -%! assert ([t(end), y(end,:)], [2, fref], 1e-3); -%!test ## anonymous function instead of real function -%! fvdb = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)]; -%! [t, y] = ode23 (fvdb, [0 2], [2 0]); -%! assert ([t(end), y(end,:)], [2, fref], 1e-3); -%!test ## extra input arguments passed through -%! [t, y] = ode23 (@fpol, [0 2], [2 0], 12, 13, "KL"); -%! assert ([t(end), y(end,:)], [2, fref], 1e-3); -%!test ## empty OdePkg structure *but* extra input arguments -%! opt = odeset; -%! [t, y] = ode23 (@fpol, [0 2], [2 0], opt, 12, 13, "KL"); -%! assert ([t(end), y(end,:)], [2, fref], 1e-2); -%!error ## strange OdePkg structure -%! opt = struct ("foo", 1); -%! [t, y] = ode23 (@fpol, [0 2], [2 0], opt); -%!test ## Solve vdp in fixed step sizes -%! opt = odeset("TimeStepSize",0.1); -%! [t, y] = ode23 (@fpol, [0,2], [2 0], opt); -%! assert (t(:), [0:0.1:2]', 1e-3); -%!test ## Solve another anonymous function below zero -%! ref = [0, 14.77810590694212]; -%! [t, y] = ode23 (@(t,y) y, [-2 0], 2); -%! assert ([t(end), y(end,:)], ref, 1e-2); -%!test ## InitialStep option -%! opt = odeset ("InitialStep", 1e-8); -%! [t, y] = ode23 (@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 = ode23 (@fpol, [0 0.2], [2 0], opt); -%! assert ([sol.x(5)-sol.x(4)], [1e-3], 1e-4); -%!test ## Solve in backward direction starting at t=0 -%! ref = [-1.205364552835178, 0.951542399860817]; -%! sol = ode23 (@fpol, [0 -2], [2 0]); -%! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 5e-3); -%!test ## Solve in backward direction starting at t=2 -%! ref = [-1.205364552835178, 0.951542399860817]; -%! sol = ode23 (@fpol, [2 0 -2], fref); -%! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 2e-2); -%!test ## Solve another anonymous function in backward direction -%! ref = [-1, 0.367879437558975]; -%! sol = ode23 (@(t,y) y, [0 -1], 1); -%! assert ([sol.x(end), sol.y(end,:)], ref, 1e-2); -%!test ## Solve another anonymous function below zero -%! ref = [0, 14.77810590694212]; -%! sol = ode23 (@(t,y) y, [-2 0], 2); -%! assert ([sol.x(end), sol.y(end,:)], ref, 1e-2); -%!test ## Solve in backward direction starting at t=0 with MaxStep option -%! ref = [-1.205364552835178, 0.951542399860817]; -%! opt = odeset ("MaxStep", 1e-3); -%! sol = ode23 (@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, ref], 1e-3); -%!test ## AbsTol option -%! opt = odeset ("AbsTol", 1e-5); -%! sol = ode23 (@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 = ode23 (@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 = ode23 (@fpol, [0 2], [2 0], opt); -%! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-4); -%!test ## Keeps initial values while integrating -%! opt = odeset ("NonNegative", 2); -%! sol = ode23 (@fpol, [0 2], [2 0], opt); -%! assert ([sol.x(end), sol.y(end,:)], [2, 2, 0], 1e-1); -%!test ## Details of OutputSel and Refine can't be tested -%! opt = odeset ("OutputFcn", @fout, "OutputSel", 1, "Refine", 5); -%! sol = ode23 (@fpol, [0 2], [2 0], opt); -%!test ## Stats must add further elements in sol -%! opt = odeset ("Stats", "on"); -%! sol = ode23 (@fpol, [0 2], [2 0], opt); -%! assert (isfield (sol, "stats")); -%! assert (isfield (sol.stats, "nsteps")); -%!test ## Events option add further elements in sol -%! opt = odeset ("Events", @feve); -%! sol = ode23 (@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 = ode23 (@fpol, [0 10], [2 0], opt); -%! assert ([sol.ie, sol.xe, sol.ye], ... -%! [2.0, 2.496110, -0.830550, -2.677589], .5e-1); -%!test ## Events option, five output arguments -%! opt = odeset ("Events", @fevn, "NormControl", "on"); -%! [t, y, vxe, ye, vie] = ode23 (@fpol, [0 10], [2 0], opt); -%! assert ([vie, vxe, ye], ... -%! [2.0, 2.496110, -0.830550, -2.677589], 1e-1); -%! -%!test ## Mass option as function -%! opt = odeset ("Mass", @fmas); -%! sol = ode23 (@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 = ode23 (@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 = ode23 (@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 = ode23 (@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 = ode23 (@fpol, [0 2], [2 0], opt); -%! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-3); +%!test # two output arguments +%! [t, y] = ode23 (@fpol, [0 2], [2 0]); +%! assert ([t(end), y(end,:)], [2, fref], 1e-3); +%!test # anonymous function instead of real function +%! fvdb = @(t,y) [y(2); (1 - y(1)^2) * y(2) - y(1)]; +%! [t, y] = ode23 (fvdb, [0 2], [2 0]); +%! assert ([t(end), y(end,:)], [2, fref], 1e-3); +%!test # extra input arguments passed through +%! [t, y] = ode23 (@fpol, [0 2], [2 0], 12, 13, "KL"); +%! assert ([t(end), y(end,:)], [2, fref], 1e-3); +%!test # empty OdePkg structure *but* extra input arguments +%! opt = odeset; +%! [t, y] = ode23 (@fpol, [0 2], [2 0], opt, 12, 13, "KL"); +%! assert ([t(end), y(end,:)], [2, fref], 1e-2); +%!test # Solve vdp in fixed step sizes +%! opt = odeset("TimeStepSize",0.1); +%! [t, y] = ode23 (@fpol, [0,2], [2 0], opt); +%! assert (t(:), [0:0.1:2]', 1e-3); +%!test # Solve another anonymous function below zero +%! ref = [0, 14.77810590694212]; +%! [t, y] = ode23 (@(t,y) y, [-2 0], 2); +%! assert ([t(end), y(end,:)], ref, 1e-2); +%!test # InitialStep option +%! opt = odeset ("InitialStep", 1e-8); +%! [t, y] = ode23 (@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 = ode23 (@fpol, [0 0.2], [2 0], opt); +%! assert ([sol.x(5)-sol.x(4)], [1e-3], 1e-4); +%!test # Solve in backward direction starting at t=0 +%! ref = [-1.205364552835178, 0.951542399860817]; +%! sol = ode23 (@fpol, [0 -2], [2 0]); +%! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 5e-3); +%!test # Solve in backward direction starting at t=2 +%! ref = [-1.205364552835178, 0.951542399860817]; +%! sol = ode23 (@fpol, [2 0 -2], fref); +%! assert ([sol.x(end), sol.y(end,:)], [-2, ref], 2e-2); +%!test # Solve another anonymous function in backward direction +%! ref = [-1, 0.367879437558975]; +%! sol = ode23 (@(t,y) y, [0 -1], 1); +%! assert ([sol.x(end), sol.y(end,:)], ref, 1e-2); +%!test # Solve another anonymous function below zero +%! ref = [0, 14.77810590694212]; +%! sol = ode23 (@(t,y) y, [-2 0], 2); +%! assert ([sol.x(end), sol.y(end,:)], ref, 1e-2); +%!test # Solve in backward direction starting at t=0 with MaxStep option +%! ref = [-1.205364552835178, 0.951542399860817]; +%! opt = odeset ("MaxStep", 1e-3); +%! sol = ode23 (@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, ref], 1e-3); +%!test # AbsTol option +%! opt = odeset ("AbsTol", 1e-5); +%! sol = ode23 (@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 = ode23 (@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 = ode23 (@fpol, [0 2], [2 0], opt); +%! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-4); +%!test # Keeps initial values while integrating +%! opt = odeset ("NonNegative", 2); +%! sol = ode23 (@fpol, [0 2], [2 0], opt); +%! assert ([sol.x(end), sol.y(end,:)], [2, 2, 0], 1e-1); +%!test # Details of OutputSel and Refine can't be tested +%! opt = odeset ("OutputFcn", @fout, "OutputSel", 1, "Refine", 5); +%! sol = ode23 (@fpol, [0 2], [2 0], opt); +%!test # Stats must add further elements in sol +%! opt = odeset ("Stats", "on"); +%! sol = ode23 (@fpol, [0 2], [2 0], opt); +%! assert (isfield (sol, "stats")); +%! assert (isfield (sol.stats, "nsteps")); +%!test # Events option add further elements in sol +%! opt = odeset ("Events", @feve); +%! sol = ode23 (@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 +%! warning ("off", "integrate_adaptive:unexpected_termination", "local"); +%! opt = odeset ("Events", @fevn, "NormControl", "on"); +%! sol = ode23 (@fpol, [0 10], [2 0], opt); +%! assert ([sol.ie, sol.xe, sol.ye], ... +%! [2.0, 2.496110, -0.830550, -2.677589], .5e-1); +%!test # Events option, five output arguments +%! warning ("off", "integrate_adaptive:unexpected_termination", "local"); +%! opt = odeset ("Events", @fevn, "NormControl", "on"); +%! [t, y, vxe, ye, vie] = ode23 (@fpol, [0 10], [2 0], opt); +%! assert ([vie, vxe, ye], [2.0, 2.496110, -0.830550, -2.677589], 1e-1); +%!test # Mass option as function +%! opt = odeset ("Mass", @fmas); +%! sol = ode23 (@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 = ode23 (@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 = ode23 (@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 = ode23 (@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 = ode23 (@fpol, [0 2], [2 0], opt); +%! assert ([sol.x(end), sol.y(end,:)], [2, fref], 1e-3); %! %! ## test for MvPattern option is missing %! ## test for InitialSlope option is missing %! ## test for MaxOrder option is missing -%! -%! warning ("on", "Octave:InvalidArgument"); -## Local Variables: *** -## mode: octave *** -## End: *** +## Test input validation +%!error ode23 () +%!error ode23 (1) +%!error ode23 (1,2) +%!error <TRANGE must be a numeric> +%! ode23 (@fpol, {[0 25]}, [3 15 1]); +%!error <TRANGE must be a .* vector> +%! ode23 (@fpol, [0 25; 25 0], [3 15 1]); +%!error <TRANGE must contain at least 2 elements> +%! ode23 (@fpol, [1], [3 15 1]); +%!error <invalid time span> +%! ode23 (@fpol, [1 1], [3 15 1]); +%!error <INIT must be a numeric> +%! ode23 (@fpol, [0 25], {[3 15 1]}); +%!error <INIT must be a .* vector> +%! ode23 (@fpol, [0 25], [3 15 1; 3 15 1]); +%!error <FUN must be a valid function handle> +%! ode23 (1, [0 25], [3 15 1]); +%!error # strange ODEOPT structure +%! opt = struct ("foo", 1); +%! [t, y] = ode23 (@fpol, [0 2], [2 0], opt); +
--- a/scripts/ode/ode45.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/ode/ode45.m Fri Mar 11 17:12:17 2016 -0800 @@ -77,7 +77,7 @@ ## @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]); +## [@var{t},@var{y}] = ode45 (fvdp, [0, 20], [2, 0]); ## @end group ## @end example ## @seealso{odeset, odeget} @@ -110,7 +110,7 @@ odeopts.funarguments = {}; endif - if (! isvector (trange) || ! isnumeric (trange)) + if (! isnumeric (trange) || ! isvector (trange)) error ("Octave:invalid-input-arg", "ode45: TRANGE must be a numeric vector"); endif @@ -129,7 +129,7 @@ endif trange = trange(:); - if (! isvector (init) || ! isnumeric (init)) + if (! isnumeric (init) || ! isvector (init)) error ("Octave:invalid-input-arg", "ode45: INIT must be a numeric vector"); endif @@ -296,7 +296,7 @@ "ode45: option 'BDF' is ignored by this solver\n"); endif - ## Starting the initialisation of the core solver ode45 + ## Starting the initialization of the core solver ode45 if (havemasshandle) # Handle only the dynamic mass matrix, if (massdependence) # constant mass matrices have already @@ -396,17 +396,17 @@ %!function ydot = fpol (t, y) # The Van der Pol %! ydot = [y(2); (1 - y(1)^2) * y(2) - y(1)]; %!endfunction -%!function ref = fref () # The computed reference solution +%!function ref = fref () # The computed reference solution %! ref = [0.32331666704577, -1.83297456798624]; %!endfunction -%!function jac = fjac (t, y, varargin) # its Jacobian +%!function jac = fjac (t, y, varargin) # its Jacobian %! jac = [0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2]; %!endfunction -%!function jac = fjcc (t, y, varargin) # sparse type +%!function jac = fjcc (t, y, varargin) # sparse type %! jac = sparse ([0, 1; -1 - 2 * y(1) * y(2), 1 - y(1)^2]) %!endfunction %!function [val, trm, dir] = feve (t, y, varargin) -%! val = fpol (t, y, varargin); # We use the derivatives +%! val = fpol (t, y, varargin); # We use the derivatives %! trm = zeros (2,1); # that's why component 2 %! dir = ones (2,1); # seems to not be exact %!endfunction @@ -434,16 +434,6 @@ %! endif %!endfunction %! -%! ## Turn off output of warning messages for all tests, -%! ## turn them on again after the last test is called -%!test -%! warning ("off", "Octave:invalid-input-arg", "local"); -%!error # first input arg is not a function -%! [t, y] = ode45 (1, [0 25], [3 15 1]); -%!error # input argument number one as name of non existing function -%! [t, y] = ode45 ("non-existing-function", [0 25], [3 15 1]); -%!error # input argument number two -%! [t, y] = ode45 (@fpol, 1, [3 15 1]); %!test # two output arguments %! [t, y] = ode45 (@fpol, [0 2], [2 0]); %! assert ([t(end), y(end,:)], [2, fref], 1e-2); @@ -464,9 +454,6 @@ %! opt = odeset; %! [t, y] = ode45 (@fpol, [0 2], [2 0], opt, 12, 13, "KL"); %! assert ([t(end), y(end,:)], [2, fref], 1e-2); -%!error # strange ODEOPT structure -%! opt = struct ("foo", 1); -%! [t, y] = ode45 (@fpol, [0 2], [2 0], opt); %!test # Solve vdp in fixed step sizes %! opt = odeset("TimeStepSize", 0.1); %! [t, y] = ode45 (@fpol, [0,2], [2 0], opt); @@ -547,11 +534,13 @@ %! assert (isfield (sol, "xe")); %! assert (isfield (sol, "ye")); %!test # Events option, now stop integration +%! warning ("off", "integrate_adaptive:unexpected_termination", "local"); %! 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], 6e-1); %!test # Events option, five output arguments +%! warning ("off", "integrate_adaptive:unexpected_termination", "local"); %! opt = odeset ("Events", @fevn, "NormControl", "on"); %! [t, y, vxe, ye, vie] = ode45 (@fpol, [0 10], [2 0], opt); %! assert ([vie, vxe, ye], @@ -576,3 +565,25 @@ %! 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); + +%!error ode45 () +%!error ode45 (1) +%!error 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 <FUN must be a valid function handle> +%! ode45 (1, [0 25], [3 15 1]); +%!error # strange ODEOPT structure +%! opt = struct ("foo", 1); +%! [t, y] = ode45 (@fpol, [0 2], [2 0], opt); +
--- a/scripts/sparse/pcg.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/sparse/pcg.m Fri Mar 11 17:12:17 2016 -0800 @@ -494,8 +494,8 @@ %!test %! ## solve small indefinite diagonal system -%! ## despite A is indefinite, the iteration continues and converges -%! ## indefiniteness of A is detected +%! ## Despite A being indefinite, the iteration continues and converges. +%! ## The indefiniteness of A is detected. %! %! N = 10; %! A = diag([1:N] .* (-ones(1, N) .^ 2)); b = rand (N, 1); @@ -519,10 +519,9 @@ %! assert (relres > 1.0); %! assert (iter, 20); # should perform max allowable default number of iterations -%!test -%! ## solve tridiagonal system with 'perfect' preconditioner -%! ## which converges in one iteration, so the eigest does not -%! ## work and issues a warning +%!warning <iteration converged too fast> +%! ## solve tridiagonal system with "perfect" preconditioner which converges +%! ## in one iteration, so the eigest does not work and issues a warning. %! %! N = 100; %! A = zeros (N, N);
--- a/scripts/special-matrix/hankel.m Sat Mar 12 00:16:18 2016 -0500 +++ b/scripts/special-matrix/hankel.m Fri Mar 11 17:12:17 2016 -0800 @@ -90,7 +90,8 @@ %!assert (hankel (1:3), [1,2,3;2,3,0;3,0,0]) %!assert (hankel (1:3,3:6), [1,2,3,4;2,3,4,5;3,4,5,6]) %!assert (hankel (1:3,3:4), [1,2;2,3;3,4]) -%!assert (hankel (1:3,4:6), [1,2,3;2,3,5;3,5,6]) +%!warning <column wins anti-diagonal conflict> +%! assert (hankel (1:3,4:6), [1,2,3;2,3,5;3,5,6]); %!error hankel () %!error hankel (1, 2, 3)