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