Mercurial > octave
annotate scripts/optimization/lsqnonneg.m @ 11588:d5bd2766c640
style fixes for warning and error messages in script files
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Thu, 20 Jan 2011 17:51:13 -0500 |
parents | c792872f8942 |
children | b9a89ca0fb75 |
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 | |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8600
diff
changeset
|
66 ## PKG_ADD: __all_opts__ ("lsqnonneg"); |
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8600
diff
changeset
|
67 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
68 ## 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
|
69 ## 161 of Solving Least Squares Problems. |
7682 | 70 |
8600 | 71 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ()) |
72 | |
73 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
74 x = optimset ("MaxIter", 1e5); | |
75 return | |
76 endif | |
77 | |
78 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options))) | |
79 print_usage (); | |
80 endif | |
7682 | 81 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
82 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 83 m = rows (c); |
84 n = columns (c); | |
7682 | 85 if (isempty (x)) |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
86 ## Initial guess is 0s. |
8600 | 87 x = zeros (n, 1); |
88 else | |
89 ## ensure nonnegative guess. | |
90 x = max (x, 0); | |
7682 | 91 endif |
92 | |
8600 | 93 useqr = m >= n; |
8507 | 94 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 95 |
8600 | 96 ## Initialize P, according to zero pattern of x. |
97 p = find (x > 0).'; | |
98 if (useqr) | |
99 ## Initialize the QR factorization, economized form. | |
100 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
101 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
102 |
7682 | 103 iter = 0; |
8600 | 104 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
105 ## LH3: test for completion. |
8600 | 106 while (iter < max_iter) |
107 while (iter < max_iter) | |
108 iter++; | |
109 | |
110 ## LH6: compute the positive matrix and find the min norm solution | |
111 ## of the positive problem. | |
112 if (useqr) | |
113 xtmp = r \ q'*d; | |
114 else | |
115 xtmp = c(:,p) \ d; | |
116 endif | |
117 idx = find (xtmp < 0); | |
118 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
119 if (isempty (idx)) |
8600 | 120 ## LH7: tmp solution found, iterate. |
121 x(:) = 0; | |
122 x(p) = xtmp; | |
123 break; | |
124 else | |
125 ## LH8, LH9: find the scaling factor. | |
126 pidx = p(idx); | |
127 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
128 alpha = min (sf); | |
129 ## LH10: adjust X. | |
130 xx = zeros (n, 1); | |
131 xx(p) = xtmp; | |
132 x += alpha*(xx - x); | |
133 ## LH11: move from P to Z all X == 0. | |
134 ## This corresponds to those indices where minimum of sf is attained. | |
135 idx = idx (sf == alpha); | |
136 p(idx) = []; | |
137 if (useqr) | |
138 ## update the QR factorization. | |
139 [q, r] = qrdelete (q, r, idx); | |
140 endif | |
141 endif | |
142 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
143 |
8600 | 144 ## compute the gradient. |
145 w = c'*(d - c*x); | |
146 w(p) = []; | |
147 if (! any (w > 0)) | |
148 if (useqr) | |
149 ## verify the solution achieved using qr updating. | |
150 ## in the best case, this should only take a single step. | |
151 useqr = false; | |
152 continue; | |
153 else | |
154 ## we're finished. | |
155 break; | |
156 endif | |
157 endif | |
158 | |
159 ## find the maximum gradient. | |
7682 | 160 idx = find (w == max (w)); |
161 if (numel (idx) > 1) | |
162 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
|
163 "a non-unique solution may be returned due to equal gradients"); |
7682 | 164 idx = idx(1); |
165 endif | |
8600 | 166 ## move the index from Z to P. Keep P sorted. |
167 z = [1:n]; z(p) = []; | |
168 zidx = z(idx); | |
169 jdx = 1 + lookup (p, zidx); | |
170 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
171 if (useqr) | |
172 ## insert the column into the QR factorization. | |
173 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
174 endif | |
7682 | 175 |
176 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
177 ## LH12: complete. |
7682 | 178 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
179 ## Generate the additional output arguments. |
7682 | 180 if (nargout > 1) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
181 resnorm = norm (c*x - d) ^ 2; |
7682 | 182 endif |
183 if (nargout > 2) | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
184 residual = d - c*x; |
7682 | 185 endif |
186 exitflag = iter; | |
8507 | 187 if (nargout > 3 && iter >= max_iter) |
7682 | 188 exitflag = 0; |
189 endif | |
190 if (nargout > 4) | |
191 output = struct ("algorithm", "nnls", "iterations", iter); | |
192 endif | |
193 if (nargout > 5) | |
8600 | 194 lambda = zeros (size (x)); |
195 lambda(p) = w; | |
7682 | 196 endif |
197 | |
198 endfunction | |
199 | |
200 ## Tests | |
201 %!test | |
202 %! C = [1 0;0 1;2 1]; | |
203 %! d = [1;3;-2]; | |
204 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps) | |
205 | |
206 %!test | |
207 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
208 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
209 %! xnew = [0;0.6929]; | |
210 %! assert (lsqnonneg (C, d), xnew, 0.0001) |