Mercurial > octave
annotate scripts/optimization/lsqnonneg.m @ 13027:b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
* fminbnd.m, fminunc.m, fsolve.m, fzero.m, lsqnonneg.m, pqpnonneg.m,
qp.m: Discard result from call to __all_opt__ in PKG_ADD command.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Mon, 29 Aug 2011 13:07:22 -0400 |
parents | d5bd2766c640 |
children | 72c96de7a403 |
rev | line source |
---|---|
11523 | 1 ## Copyright (C) 2008-2011 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}) |
24 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}] =} lsqnonneg (@dots{}) | |
25 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}] =} lsqnonneg (@dots{}) | |
26 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}] =} lsqnonneg (@dots{}) | |
27 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}] =} lsqnonneg (@dots{}) | |
28 ## @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
|
29 ## 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
|
30 ## @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
|
31 ## optional initial guess for @var{x}. |
7682 | 32 ## |
33 ## Outputs: | |
34 ## @itemize @bullet | |
35 ## @item resnorm | |
36 ## | |
37 ## 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
|
38 ## |
7682 | 39 ## @item residual |
40 ## | |
41 ## The residual: @var{d}-@var{c}*@var{x} | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
42 ## |
7682 | 43 ## @item exitflag |
44 ## | |
45 ## An indicator of convergence. 0 indicates that the iteration count | |
46 ## was exceeded, and therefore convergence was not reached; >0 indicates | |
47 ## that the algorithm converged. (The algorithm is stable and will | |
48 ## converge given enough iterations.) | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
49 ## |
7682 | 50 ## @item output |
51 ## | |
52 ## A structure with two fields: | |
53 ## @itemize @bullet | |
54 ## @item "algorithm": The algorithm used ("nnls") | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## |
7682 | 56 ## @item "iterations": The number of iterations taken. |
57 ## @end itemize | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
58 ## |
7682 | 59 ## @item lambda |
60 ## | |
61 ## Not implemented. | |
62 ## @end itemize | |
9635 | 63 ## @seealso{optimset, pqpnonneg} |
7682 | 64 ## @end deftypefn |
65 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
66 ## 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
|
67 ## 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
|
68 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
69 ## 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
|
70 ## 161 of Solving Least Squares Problems. |
7682 | 71 |
8600 | 72 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ()) |
73 | |
74 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
75 x = optimset ("MaxIter", 1e5); | |
76 return | |
77 endif | |
78 | |
79 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options))) | |
80 print_usage (); | |
81 endif | |
7682 | 82 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
83 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 84 m = rows (c); |
85 n = columns (c); | |
7682 | 86 if (isempty (x)) |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
87 ## Initial guess is 0s. |
8600 | 88 x = zeros (n, 1); |
89 else | |
90 ## ensure nonnegative guess. | |
91 x = max (x, 0); | |
7682 | 92 endif |
93 | |
8600 | 94 useqr = m >= n; |
8507 | 95 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 96 |
8600 | 97 ## Initialize P, according to zero pattern of x. |
98 p = find (x > 0).'; | |
99 if (useqr) | |
100 ## Initialize the QR factorization, economized form. | |
101 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
102 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
103 |
7682 | 104 iter = 0; |
8600 | 105 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
106 ## LH3: test for completion. |
8600 | 107 while (iter < max_iter) |
108 while (iter < max_iter) | |
109 iter++; | |
110 | |
111 ## LH6: compute the positive matrix and find the min norm solution | |
112 ## of the positive problem. | |
113 if (useqr) | |
114 xtmp = r \ q'*d; | |
115 else | |
116 xtmp = c(:,p) \ d; | |
117 endif | |
118 idx = find (xtmp < 0); | |
119 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
120 if (isempty (idx)) |
8600 | 121 ## LH7: tmp solution found, iterate. |
122 x(:) = 0; | |
123 x(p) = xtmp; | |
124 break; | |
125 else | |
126 ## LH8, LH9: find the scaling factor. | |
127 pidx = p(idx); | |
128 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
129 alpha = min (sf); | |
130 ## LH10: adjust X. | |
131 xx = zeros (n, 1); | |
132 xx(p) = xtmp; | |
133 x += alpha*(xx - x); | |
134 ## LH11: move from P to Z all X == 0. | |
135 ## This corresponds to those indices where minimum of sf is attained. | |
136 idx = idx (sf == alpha); | |
137 p(idx) = []; | |
138 if (useqr) | |
139 ## update the QR factorization. | |
140 [q, r] = qrdelete (q, r, idx); | |
141 endif | |
142 endif | |
143 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
144 |
8600 | 145 ## compute the gradient. |
146 w = c'*(d - c*x); | |
147 w(p) = []; | |
148 if (! any (w > 0)) | |
149 if (useqr) | |
150 ## verify the solution achieved using qr updating. | |
151 ## in the best case, this should only take a single step. | |
152 useqr = false; | |
153 continue; | |
154 else | |
155 ## we're finished. | |
156 break; | |
157 endif | |
158 endif | |
159 | |
160 ## find the maximum gradient. | |
7682 | 161 idx = find (w == max (w)); |
162 if (numel (idx) > 1) | |
163 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
|
164 "a non-unique solution may be returned due to equal gradients"); |
7682 | 165 idx = idx(1); |
166 endif | |
8600 | 167 ## move the index from Z to P. Keep P sorted. |
168 z = [1:n]; z(p) = []; | |
169 zidx = z(idx); | |
170 jdx = 1 + lookup (p, zidx); | |
171 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
172 if (useqr) | |
173 ## insert the column into the QR factorization. | |
174 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
175 endif | |
7682 | 176 |
177 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
178 ## LH12: complete. |
7682 | 179 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
180 ## Generate the additional output arguments. |
7682 | 181 if (nargout > 1) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
182 resnorm = norm (c*x - d) ^ 2; |
7682 | 183 endif |
184 if (nargout > 2) | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
185 residual = d - c*x; |
7682 | 186 endif |
187 exitflag = iter; | |
8507 | 188 if (nargout > 3 && iter >= max_iter) |
7682 | 189 exitflag = 0; |
190 endif | |
191 if (nargout > 4) | |
192 output = struct ("algorithm", "nnls", "iterations", iter); | |
193 endif | |
194 if (nargout > 5) | |
8600 | 195 lambda = zeros (size (x)); |
196 lambda(p) = w; | |
7682 | 197 endif |
198 | |
199 endfunction | |
200 | |
201 ## Tests | |
202 %!test | |
203 %! C = [1 0;0 1;2 1]; | |
204 %! d = [1;3;-2]; | |
205 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps) | |
206 | |
207 %!test | |
208 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
209 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
210 %! xnew = [0;0.6929]; | |
211 %! assert (lsqnonneg (C, d), xnew, 0.0001) |