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