Mercurial > octave-nkf
annotate scripts/optimization/lsqnonneg.m @ 19867:9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Try to trim long lines to < 80 chars.
Use '##' for single line comments.
Use '(...)' around tests for if/elseif/switch/while.
Abut cell indexing operator '{' next to variable.
Abut array indexing operator '(' next to variable.
Use space between negation operator '!' and following expression.
Use two newlines between endfunction and start of %!test or %!demo code.
Remove unnecessary parens grouping between short-circuit operators.
Remove stray extra spaces (typos) between variables and assignment operators.
Remove stray extra spaces from ends of lines.
author | Rik <rik@octave.org> |
---|---|
date | Mon, 23 Feb 2015 14:54:39 -0800 |
parents | 4197fc428c7d |
children | f1d0f506ee78 |
rev | line source |
---|---|
19731
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19105
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 |
10846
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
31 ## @code{@var{x} >= 0}. @var{c} and @var{d} must be real. @var{x0} is an |
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
32 ## optional initial guess for @var{x}. |
14787
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
33 ## Currently, @code{lsqnonneg} |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
34 ## recognizes these options: @qcode{"MaxIter"}, @qcode{"TolX"}. |
17097
e7a059a9a644
doc: Use XREF as anchor prefix in documentation for clearer results in Info viewer.
Rik <rik@octave.org>
parents:
16772
diff
changeset
|
35 ## For a description of these options, see @ref{XREFoptimset,,optimset}. |
7682 | 36 ## |
37 ## Outputs: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
38 ## |
7682 | 39 ## @itemize @bullet |
40 ## @item resnorm | |
41 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14787
diff
changeset
|
42 ## 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
|
43 ## |
7682 | 44 ## @item residual |
45 ## | |
46 ## The residual: @var{d}-@var{c}*@var{x} | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
47 ## |
7682 | 48 ## @item exitflag |
49 ## | |
50 ## An indicator of convergence. 0 indicates that the iteration count | |
51 ## was exceeded, and therefore convergence was not reached; >0 indicates | |
52 ## that the algorithm converged. (The algorithm is stable and will | |
53 ## converge given enough iterations.) | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
54 ## |
7682 | 55 ## @item output |
56 ## | |
57 ## A structure with two fields: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
58 ## |
7682 | 59 ## @itemize @bullet |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
60 ## @item @qcode{"algorithm"}: The algorithm used (@qcode{"nnls"}) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
61 ## |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17097
diff
changeset
|
62 ## @item @qcode{"iterations"}: The number of iterations taken. |
7682 | 63 ## @end itemize |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
64 ## |
7682 | 65 ## @item lambda |
66 ## | |
67 ## Not implemented. | |
68 ## @end itemize | |
19105
c9113e28fae8
New lscov function (bug #43118)
Pantxo Diribarne <pantxo.diribarne@gmail.com>
parents:
17744
diff
changeset
|
69 ## @seealso{optimset, pqpnonneg, lscov} |
7682 | 70 ## @end deftypefn |
71 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
72 ## 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
|
73 ## 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
|
74 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
75 ## 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
|
76 ## 161 of Solving Least Squares Problems. |
7682 | 77 |
8600 | 78 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ()) |
79 | |
80 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
81 x = optimset ("MaxIter", 1e5); | |
17312
088d014a7fe2
Use semicolon after "return" statement in core m-files.
Rik <rik@octave.org>
parents:
17281
diff
changeset
|
82 return; |
8600 | 83 endif |
84 | |
19867
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
85 if (nargin < 2 || nargin > 4 |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
86 || ! (ismatrix (c) && ismatrix (d) && isstruct (options))) |
8600 | 87 print_usage (); |
88 endif | |
7682 | 89 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
90 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 91 m = rows (c); |
92 n = columns (c); | |
7682 | 93 if (isempty (x)) |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
94 ## Initial guess is 0s. |
8600 | 95 x = zeros (n, 1); |
96 else | |
97 ## ensure nonnegative guess. | |
98 x = max (x, 0); | |
7682 | 99 endif |
100 | |
8600 | 101 useqr = m >= n; |
8507 | 102 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 103 |
8600 | 104 ## Initialize P, according to zero pattern of x. |
105 p = find (x > 0).'; | |
106 if (useqr) | |
107 ## Initialize the QR factorization, economized form. | |
108 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
109 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
110 |
7682 | 111 iter = 0; |
8600 | 112 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
113 ## LH3: test for completion. |
8600 | 114 while (iter < max_iter) |
115 while (iter < max_iter) | |
116 iter++; | |
117 | |
118 ## LH6: compute the positive matrix and find the min norm solution | |
119 ## of the positive problem. | |
120 if (useqr) | |
121 xtmp = r \ q'*d; | |
122 else | |
123 xtmp = c(:,p) \ d; | |
124 endif | |
125 idx = find (xtmp < 0); | |
126 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
127 if (isempty (idx)) |
8600 | 128 ## LH7: tmp solution found, iterate. |
129 x(:) = 0; | |
130 x(p) = xtmp; | |
131 break; | |
132 else | |
133 ## LH8, LH9: find the scaling factor. | |
134 pidx = p(idx); | |
135 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
136 alpha = min (sf); | |
137 ## LH10: adjust X. | |
138 xx = zeros (n, 1); | |
139 xx(p) = xtmp; | |
140 x += alpha*(xx - x); | |
141 ## LH11: move from P to Z all X == 0. | |
142 ## This corresponds to those indices where minimum of sf is attained. | |
143 idx = idx (sf == alpha); | |
144 p(idx) = []; | |
145 if (useqr) | |
146 ## update the QR factorization. | |
147 [q, r] = qrdelete (q, r, idx); | |
148 endif | |
149 endif | |
150 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
151 |
8600 | 152 ## compute the gradient. |
153 w = c'*(d - c*x); | |
154 w(p) = []; | |
14787
acb09716fc94
lsqnonneg have tolerance option for convergence (bug #33347)
Axel Mathéi <axel.mathei@gmail.com>
parents:
14366
diff
changeset
|
155 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
|
156 if (! any (w > tolx)) |
8600 | 157 if (useqr) |
158 ## verify the solution achieved using qr updating. | |
159 ## in the best case, this should only take a single step. | |
160 useqr = false; | |
161 continue; | |
162 else | |
163 ## we're finished. | |
164 break; | |
165 endif | |
166 endif | |
167 | |
168 ## find the maximum gradient. | |
7682 | 169 idx = find (w == max (w)); |
170 if (numel (idx) > 1) | |
171 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
|
172 "a non-unique solution may be returned due to equal gradients"); |
7682 | 173 idx = idx(1); |
174 endif | |
8600 | 175 ## move the index from Z to P. Keep P sorted. |
176 z = [1:n]; z(p) = []; | |
177 zidx = z(idx); | |
178 jdx = 1 + lookup (p, zidx); | |
179 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
180 if (useqr) | |
181 ## insert the column into the QR factorization. | |
182 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
183 endif | |
7682 | 184 |
185 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
186 ## LH12: complete. |
7682 | 187 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
188 ## Generate the additional output arguments. |
7682 | 189 if (nargout > 1) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
190 resnorm = norm (c*x - d) ^ 2; |
7682 | 191 endif |
192 if (nargout > 2) | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
193 residual = d - c*x; |
7682 | 194 endif |
195 exitflag = iter; | |
8507 | 196 if (nargout > 3 && iter >= max_iter) |
7682 | 197 exitflag = 0; |
198 endif | |
199 if (nargout > 4) | |
200 output = struct ("algorithm", "nnls", "iterations", iter); | |
201 endif | |
202 if (nargout > 5) | |
8600 | 203 lambda = zeros (size (x)); |
204 lambda(p) = w; | |
7682 | 205 endif |
206 | |
207 endfunction | |
208 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
209 |
7682 | 210 %!test |
211 %! C = [1 0;0 1;2 1]; | |
212 %! 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
|
213 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps); |
7682 | 214 |
215 %!test | |
216 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
217 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
218 %! 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
|
219 %! 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
|
220 |