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