Mercurial > octave
annotate scripts/optimization/lsqnonneg.m @ 21171:2935d56203a4 stable
Fix regressions caused by ismatrix definition change (partial fix bug #47036).
* inputdlg.m: Test that linespec isnumeric.
* uigetfile.m: Check that position property value isnumeric.
* fminunc.m: Check that x0 isnumeric.
* fsolve.m: Check that x0 isnumeric.
* lsqnonneg.m: Check that inputs C & D are both isnumeric and ismatrix.
* pqpnonneg.m: Check that inputs C & D are both isnumeric and ismatrix.
* bicg.m: Check input A issquare. Rephrase error messages.
* bicgstab.m: Check input A issquare. Rephrase error messages.
* cgs.m: Check input A issquare. Rephrase error messages.
* gmres.m: Check input A issquare. Rephrase error messages.
Change BIST test to match new error message.
* qmr.m: Check input A issquare. Rephrase error messages.
* spconvert.m: Check nargin first. Simplify input validation.
Wrap long error message to < 80 chars.
* treeplot.m: Simplify input validation.
author | Rik <rik@octave.org> |
---|---|
date | Mon, 01 Feb 2016 22:59:43 -0800 |
parents | f1d0f506ee78 |
children | 3be6a07e8bad |
rev | line source |
---|---|
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19071
diff
changeset
|
1 ## Copyright (C) 2008-2015 Bill Denney |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
2 ## Copyright (C) 2008 Jaroslav Hajek |
8600 | 3 ## Copyright (C) 2009 VZLU Prague |
7682 | 4 ## |
5 ## This file is part of Octave. | |
6 ## | |
7 ## Octave is free software; you can redistribute it and/or modify it | |
8 ## under the terms of the GNU General Public License as published by | |
9 ## the Free Software Foundation; either version 3 of the License, or (at | |
10 ## your option) any later version. | |
11 ## | |
12 ## Octave is distributed in the hope that it will be useful, but | |
13 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
14 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU | |
15 ## General Public License for more details. | |
16 ## | |
17 ## You should have received a copy of the GNU General Public License | |
18 ## along with Octave; see the file COPYING. If not, see | |
19 ## <http://www.gnu.org/licenses/>. | |
20 | |
21 ## -*- texinfo -*- | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
9635
diff
changeset
|
22 ## @deftypefn {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}) |
7682 | 23 ## @deftypefnx {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0}) |
14787
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
24 ## @deftypefnx {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0}, @var{options}) |
7682 | 25 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}] =} lsqnonneg (@dots{}) |
26 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}] =} lsqnonneg (@dots{}) | |
27 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}] =} lsqnonneg (@dots{}) | |
28 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}] =} lsqnonneg (@dots{}) | |
29 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}, @var{lambda}] =} lsqnonneg (@dots{}) | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
30 ## Minimize @code{norm (@var{c}*@var{x} - d)} subject to |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
31 ## @code{@var{x} >= 0}. |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
32 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
33 ## @var{c} and @var{d} must be real. |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
34 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
35 ## @var{x0} is an optional initial guess for @var{x}. |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
36 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
37 ## Currently, @code{lsqnonneg} recognizes these options: @qcode{"MaxIter"}, |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
38 ## @qcode{"TolX"}. For a description of these options, see |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
39 ## @ref{XREFoptimset,,optimset}. |
7682 | 40 ## |
41 ## Outputs: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
42 ## |
7682 | 43 ## @itemize @bullet |
44 ## @item resnorm | |
45 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14787
diff
changeset
|
46 ## The squared 2-norm of the residual: norm (@var{c}*@var{x}-@var{d})^2 |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
47 ## |
7682 | 48 ## @item residual |
49 ## | |
50 ## The residual: @var{d}-@var{c}*@var{x} | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
51 ## |
7682 | 52 ## @item exitflag |
53 ## | |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
54 ## An indicator of convergence. 0 indicates that the iteration count was |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
55 ## exceeded, and therefore convergence was not reached; >0 indicates that the |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
56 ## algorithm converged. (The algorithm is stable and will converge given |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
57 ## enough iterations.) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
58 ## |
7682 | 59 ## @item output |
60 ## | |
61 ## A structure with two fields: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
62 ## |
7682 | 63 ## @itemize @bullet |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
64 ## @item @qcode{"algorithm"}: The algorithm used (@qcode{"nnls"}) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
65 ## |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
66 ## @item @qcode{"iterations"}: The number of iterations taken. |
7682 | 67 ## @end itemize |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
68 ## |
7682 | 69 ## @item lambda |
70 ## | |
71 ## Not implemented. | |
72 ## @end itemize | |
19071
c9113e28fae8
New lscov function (bug #43118)
Pantxo Diribarne <pantxo.diribarne@gmail.com>
parents:
17744
diff
changeset
|
73 ## @seealso{optimset, pqpnonneg, lscov} |
7682 | 74 ## @end deftypefn |
75 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
76 ## PKG_ADD: ## Discard result to avoid polluting workspace with ans at startup. |
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
77 ## PKG_ADD: [~] = __all_opts__ ("lsqnonneg"); |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8600
diff
changeset
|
78 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
79 ## This is implemented from Lawson and Hanson's 1973 algorithm on page |
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
80 ## 161 of Solving Least Squares Problems. |
7682 | 81 |
8600 | 82 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ()) |
83 | |
84 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
85 x = optimset ("MaxIter", 1e5); | |
17312
088d014a7fe2
Use semicolon after "return" statement in core m-files.
Rik <rik@octave.org>
parents:
17281
diff
changeset
|
86 return; |
8600 | 87 endif |
88 | |
19833
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
89 if (nargin < 2 || nargin > 4 |
21171
2935d56203a4
Fix regressions caused by ismatrix definition change (partial fix bug #47036).
Rik <rik@octave.org>
parents:
20165
diff
changeset
|
90 || ! (isnumeric (c) && ismatrix (c)) |
2935d56203a4
Fix regressions caused by ismatrix definition change (partial fix bug #47036).
Rik <rik@octave.org>
parents:
20165
diff
changeset
|
91 || ! (isnumeric (d) && ismatrix (d)) |
2935d56203a4
Fix regressions caused by ismatrix definition change (partial fix bug #47036).
Rik <rik@octave.org>
parents:
20165
diff
changeset
|
92 || ! isstruct (options)) |
8600 | 93 print_usage (); |
94 endif | |
7682 | 95 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
96 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 97 m = rows (c); |
98 n = columns (c); | |
7682 | 99 if (isempty (x)) |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
100 ## Initial guess is 0s. |
8600 | 101 x = zeros (n, 1); |
102 else | |
103 ## ensure nonnegative guess. | |
104 x = max (x, 0); | |
7682 | 105 endif |
106 | |
8600 | 107 useqr = m >= n; |
8507 | 108 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 109 |
8600 | 110 ## Initialize P, according to zero pattern of x. |
111 p = find (x > 0).'; | |
112 if (useqr) | |
113 ## Initialize the QR factorization, economized form. | |
114 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
115 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
116 |
7682 | 117 iter = 0; |
8600 | 118 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
119 ## LH3: test for completion. |
8600 | 120 while (iter < max_iter) |
121 while (iter < max_iter) | |
122 iter++; | |
123 | |
124 ## LH6: compute the positive matrix and find the min norm solution | |
125 ## of the positive problem. | |
126 if (useqr) | |
127 xtmp = r \ q'*d; | |
128 else | |
129 xtmp = c(:,p) \ d; | |
130 endif | |
131 idx = find (xtmp < 0); | |
132 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
133 if (isempty (idx)) |
8600 | 134 ## LH7: tmp solution found, iterate. |
135 x(:) = 0; | |
136 x(p) = xtmp; | |
137 break; | |
138 else | |
139 ## LH8, LH9: find the scaling factor. | |
140 pidx = p(idx); | |
141 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
142 alpha = min (sf); | |
143 ## LH10: adjust X. | |
144 xx = zeros (n, 1); | |
145 xx(p) = xtmp; | |
146 x += alpha*(xx - x); | |
147 ## LH11: move from P to Z all X == 0. | |
148 ## This corresponds to those indices where minimum of sf is attained. | |
149 idx = idx (sf == alpha); | |
150 p(idx) = []; | |
151 if (useqr) | |
152 ## update the QR factorization. | |
153 [q, r] = qrdelete (q, r, idx); | |
154 endif | |
155 endif | |
156 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
157 |
8600 | 158 ## compute the gradient. |
159 w = c'*(d - c*x); | |
160 w(p) = []; | |
14787
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
161 tolx = optimget (options, "TolX", 10*eps*norm (c, 1)*length (c)); |
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
162 if (! any (w > tolx)) |
8600 | 163 if (useqr) |
164 ## verify the solution achieved using qr updating. | |
165 ## in the best case, this should only take a single step. | |
166 useqr = false; | |
167 continue; | |
168 else | |
169 ## we're finished. | |
170 break; | |
171 endif | |
172 endif | |
173 | |
174 ## find the maximum gradient. | |
7682 | 175 idx = find (w == max (w)); |
176 if (numel (idx) > 1) | |
177 warning ("lsqnonneg:nonunique", | |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
178 "a non-unique solution may be returned due to equal gradients"); |
7682 | 179 idx = idx(1); |
180 endif | |
8600 | 181 ## move the index from Z to P. Keep P sorted. |
182 z = [1:n]; z(p) = []; | |
183 zidx = z(idx); | |
184 jdx = 1 + lookup (p, zidx); | |
185 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
186 if (useqr) | |
187 ## insert the column into the QR factorization. | |
188 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
189 endif | |
7682 | 190 |
191 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
192 ## LH12: complete. |
7682 | 193 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
194 ## Generate the additional output arguments. |
7682 | 195 if (nargout > 1) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
196 resnorm = norm (c*x - d) ^ 2; |
7682 | 197 endif |
198 if (nargout > 2) | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
199 residual = d - c*x; |
7682 | 200 endif |
201 exitflag = iter; | |
8507 | 202 if (nargout > 3 && iter >= max_iter) |
7682 | 203 exitflag = 0; |
204 endif | |
205 if (nargout > 4) | |
206 output = struct ("algorithm", "nnls", "iterations", iter); | |
207 endif | |
208 if (nargout > 5) | |
8600 | 209 lambda = zeros (size (x)); |
210 lambda(p) = w; | |
7682 | 211 endif |
212 | |
213 endfunction | |
214 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
215 |
7682 | 216 %!test |
217 %! C = [1 0;0 1;2 1]; | |
218 %! d = [1;3;-2]; | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
219 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps); |
7682 | 220 |
221 %!test | |
222 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
223 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
224 %! xnew = [0;0.6929]; | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
225 %! assert (lsqnonneg (C, d), xnew, 0.0001); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
226 |