Mercurial > octave
annotate scripts/optimization/lsqnonneg.m @ 29359:7854d5752dd2
maint: merge stable to default.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 10 Feb 2021 10:10:40 -0500 |
parents | 5394d688d456 0a5b15007766 |
children | 01de0045b2e3 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
29358
0a5b15007766
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
27923
diff
changeset
|
3 ## Copyright (C) 2008-2021 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
7682 | 7 ## |
8 ## This file is part of Octave. | |
9 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
7682 | 11 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
13 ## (at your option) any later version. |
7682 | 14 ## |
15 ## Octave is distributed in the hope that it will be useful, but | |
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
18 ## GNU General Public License for more details. |
7682 | 19 ## |
20 ## You should have received a copy of the GNU General Public License | |
21 ## along with Octave; see the file COPYING. If not, see | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
7682 | 25 |
26 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
27 ## @deftypefn {} {@var{x} =} lsqnonneg (@var{c}, @var{d}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
28 ## @deftypefnx {} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
29 ## @deftypefnx {} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0}, @var{options}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
30 ## @deftypefnx {} {[@var{x}, @var{resnorm}] =} lsqnonneg (@dots{}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
31 ## @deftypefnx {} {[@var{x}, @var{resnorm}, @var{residual}] =} lsqnonneg (@dots{}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
32 ## @deftypefnx {} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}] =} lsqnonneg (@dots{}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
33 ## @deftypefnx {} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}] =} lsqnonneg (@dots{}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20735
diff
changeset
|
34 ## @deftypefnx {} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}, @var{lambda}] =} lsqnonneg (@dots{}) |
22769 | 35 ## |
36 ## Minimize @code{norm (@var{c}*@var{x} - @var{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
|
37 ## @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
|
38 ## |
22769 | 39 ## @var{c} and @var{d} must be real matrices. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
40 ## |
22769 | 41 ## @var{x0} is an optional initial guess for the solution @var{x}. |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
42 ## |
22769 | 43 ## @var{options} is an options structure to change the behavior of the |
28959
5394d688d456
doc: Use @code{} within alternate text for @xref,@pxref macros for better Info display.
Rik <rik@octave.org>
parents:
28942
diff
changeset
|
44 ## algorithm (@pxref{XREFoptimset,,@code{optimset}}). @code{lsqnonneg} |
5394d688d456
doc: Use @code{} within alternate text for @xref,@pxref macros for better Info display.
Rik <rik@octave.org>
parents:
28942
diff
changeset
|
45 ## recognizes these options: @qcode{"MaxIter"}, @qcode{"TolX"}. |
7682 | 46 ## |
47 ## Outputs: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
48 ## |
22770
31b86d30874d
lsqnonneg.m: Change @itemize to @table so docs compile.
Rik <rik@octave.org>
parents:
22769
diff
changeset
|
49 ## @table @var |
7682 | 50 ## @item resnorm |
22769 | 51 ## The squared 2-norm of the residual: @code{norm (@var{c}*@var{x}-@var{d})^2} |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
52 ## |
7682 | 53 ## @item residual |
22769 | 54 ## The residual: @code{@var{d}-@var{c}*@var{x}} |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## |
7682 | 56 ## @item exitflag |
20165
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19833
diff
changeset
|
57 ## 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
|
58 ## 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
|
59 ## 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
|
60 ## enough iterations.) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
61 ## |
7682 | 62 ## @item output |
63 ## A structure with two fields: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
64 ## |
7682 | 65 ## @itemize @bullet |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
66 ## @item @qcode{"algorithm"}: The algorithm used (@qcode{"nnls"}) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
67 ## |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
68 ## @item @qcode{"iterations"}: The number of iterations taken. |
7682 | 69 ## @end itemize |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
70 ## |
7682 | 71 ## @item lambda |
22769 | 72 ## @c FIXME: Something is output from the function, but what is it? |
73 ## Undocumented output | |
74 ## @end table | |
75 ## @seealso{pqpnonneg, lscov, optimset} | |
7682 | 76 ## @end deftypefn |
77 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
78 ## 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
|
79 ## 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
|
80 |
22769 | 81 ## This is implemented from Lawson and Hanson's 1973 algorithm on page 161 of |
82 ## Solving Least Squares Problems. | |
7682 | 83 |
22769 | 84 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, |
85 x0 = [], | |
86 options = struct ()) | |
8600 | 87 |
22769 | 88 ## Special case: called to find default optimization options |
89 if (nargin == 1 && ischar (c) && strcmp (c, "defaults")) | |
26138
804e18e3e320
Reenable query of optimization options (bugs #54952 and #55089).
Kai T. Ohlhus <k.ohlhus@gmail.com>
parents:
25579
diff
changeset
|
90 x = struct ("MaxIter", 1e5); |
17312
088d014a7fe2
Use semicolon after "return" statement in core m-files.
Rik <rik@octave.org>
parents:
17281
diff
changeset
|
91 return; |
8600 | 92 endif |
93 | |
28789
28de41192f3c
Eliminate unneeded verification of nargin, nargout in m-files.
Rik <rik@octave.org>
parents:
27923
diff
changeset
|
94 if (nargin < 2) |
8600 | 95 print_usage (); |
96 endif | |
7682 | 97 |
22769 | 98 if (! (isnumeric (c) && ismatrix (c)) || ! (isnumeric (d) && ismatrix (d))) |
99 error ("lsqnonneg: C and D must be numeric matrices"); | |
100 endif | |
101 | |
102 if (! isstruct (options)) | |
103 error ("lsqnonneg: OPTIONS must be a struct"); | |
104 endif | |
105 | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
106 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 107 m = rows (c); |
108 n = columns (c); | |
22769 | 109 if (isempty (x0)) |
110 ## Initial guess is all zeros. | |
8600 | 111 x = zeros (n, 1); |
112 else | |
113 ## ensure nonnegative guess. | |
22769 | 114 x = max (x0, 0); |
7682 | 115 endif |
116 | |
22769 | 117 useqr = (m >= n); |
8507 | 118 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 119 |
8600 | 120 ## Initialize P, according to zero pattern of x. |
121 p = find (x > 0).'; | |
122 if (useqr) | |
123 ## Initialize the QR factorization, economized form. | |
124 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
125 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
126 |
7682 | 127 iter = 0; |
8600 | 128 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
129 ## LH3: test for completion. |
8600 | 130 while (iter < max_iter) |
131 while (iter < max_iter) | |
20735
418ae0cb752f
Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents:
20231
diff
changeset
|
132 iter += 1; |
8600 | 133 |
134 ## LH6: compute the positive matrix and find the min norm solution | |
135 ## of the positive problem. | |
136 if (useqr) | |
137 xtmp = r \ q'*d; | |
138 else | |
139 xtmp = c(:,p) \ d; | |
140 endif | |
141 idx = find (xtmp < 0); | |
142 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
143 if (isempty (idx)) |
8600 | 144 ## LH7: tmp solution found, iterate. |
145 x(:) = 0; | |
146 x(p) = xtmp; | |
147 break; | |
148 else | |
149 ## LH8, LH9: find the scaling factor. | |
150 pidx = p(idx); | |
22769 | 151 sf = x(pidx) ./ (x(pidx) - xtmp(idx)); |
8600 | 152 alpha = min (sf); |
153 ## LH10: adjust X. | |
154 xx = zeros (n, 1); | |
155 xx(p) = xtmp; | |
156 x += alpha*(xx - x); | |
157 ## LH11: move from P to Z all X == 0. | |
158 ## This corresponds to those indices where minimum of sf is attained. | |
20231
83792dd9bcc1
Use in-place operators in m-files where possible.
Rik <rik@octave.org>
parents:
20165
diff
changeset
|
159 idx = idx(sf == alpha); |
8600 | 160 p(idx) = []; |
161 if (useqr) | |
162 ## update the QR factorization. | |
163 [q, r] = qrdelete (q, r, idx); | |
164 endif | |
165 endif | |
166 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
167 |
8600 | 168 ## compute the gradient. |
169 w = c'*(d - c*x); | |
170 w(p) = []; | |
14787
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
171 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
|
172 if (! any (w > tolx)) |
8600 | 173 if (useqr) |
174 ## verify the solution achieved using qr updating. | |
175 ## in the best case, this should only take a single step. | |
176 useqr = false; | |
177 continue; | |
178 else | |
179 ## we're finished. | |
180 break; | |
181 endif | |
182 endif | |
183 | |
184 ## find the maximum gradient. | |
7682 | 185 idx = find (w == max (w)); |
186 if (numel (idx) > 1) | |
187 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
|
188 "a non-unique solution may be returned due to equal gradients"); |
7682 | 189 idx = idx(1); |
190 endif | |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21178
diff
changeset
|
191 ## move the index from Z to P. Keep P sorted. |
8600 | 192 z = [1:n]; z(p) = []; |
193 zidx = z(idx); | |
194 jdx = 1 + lookup (p, zidx); | |
195 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
196 if (useqr) | |
197 ## insert the column into the QR factorization. | |
198 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
199 endif | |
7682 | 200 |
201 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
202 ## LH12: complete. |
7682 | 203 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
204 ## Generate the additional output arguments. |
22769 | 205 if (isargout (2)) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
206 resnorm = norm (c*x - d) ^ 2; |
7682 | 207 endif |
22769 | 208 if (isargout (3)) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
209 residual = d - c*x; |
7682 | 210 endif |
22769 | 211 if (isargout (4)) |
212 if (iter >= max_iter) | |
213 exitflag = 0; | |
214 else | |
215 exitflag = iter; | |
216 endif | |
7682 | 217 endif |
22769 | 218 if (isargout (5)) |
7682 | 219 output = struct ("algorithm", "nnls", "iterations", iter); |
220 endif | |
22769 | 221 if (isargout (6)) |
8600 | 222 lambda = zeros (size (x)); |
223 lambda(p) = w; | |
7682 | 224 endif |
225 | |
226 endfunction | |
227 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
228 |
7682 | 229 %!test |
230 %! C = [1 0;0 1;2 1]; | |
231 %! 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
|
232 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps); |
7682 | 233 |
234 %!test | |
235 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
236 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
237 %! 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
|
238 %! 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
|
239 |
28942
fc4bb4bd1d5e
maint: Use '##' as lead-in for full-line comments.
Rik <rik@octave.org>
parents:
28896
diff
changeset
|
240 ## Test input validation |
28896
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
241 %!error <Invalid call> lsqnonneg () |
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
242 %!error <Invalid call> lsqnonneg (1) |
22769 | 243 %!error <C .* must be numeric matrices> lsqnonneg ({1},2) |
244 %!error <C .* must be numeric matrices> lsqnonneg (ones (2,2,2),2) | |
245 %!error <D must be numeric matrices> lsqnonneg (1,{2}) | |
246 %!error <D must be numeric matrices> lsqnonneg (1, ones (2,2,2)) | |
247 %!error <OPTIONS must be a struct> lsqnonneg (1, 2, [], 3) |