changelog shortlog tags changeset files revisions annotate raw

scripts/optimization/lsqnonneg.m

changeset 10289: 4b124317dc38
parent:1bf0ce0930be
author: John W. Eaton <jwe@octave.org>
date: Tue Feb 09 20:58:55 2010 -0500 (69 minutes ago)
permissions: -rw-r--r--
description: base_properties::set_children: account for hidden children
1## Copyright (C) 2008 Bill Denney
2## Copyright (C) 2008 Jaroslav Hajek
3## Copyright (C) 2009 VZLU Prague
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, pqpnonneg}
59## @end deftypefn
60
61## PKG_ADD: __all_opts__ ("lsqnonneg");
62
63## This is implemented from Lawson and Hanson's 1973 algorithm on page
64## 161 of Solving Least Squares Problems.
65
66function [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
76
77 ## Lawson-Hanson Step 1 (LH1): initialize the variables.
78 m = rows (c);
79 n = columns (c);
80 if (isempty (x))
81 ## Initial guess is 0s.
82 x = zeros (n, 1);
83 else
84 ## ensure nonnegative guess.
85 x = max (x, 0);
86 endif
87
88 useqr = m >= n;
89 max_iter = optimget (options, "MaxIter", 1e5);
90
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);
96 endif
97
98 iter = 0;
99
100 ## LH3: test for completion.
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.
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
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
170
171 endwhile
172 ## LH12: complete.
173
174 ## Generate the additional output arguments.
175 if (nargout > 1)
176 resnorm = norm (c*x - d) ^ 2;
177 endif
178 if (nargout > 2)
179 residual = d - c*x;
180 endif
181 exitflag = iter;
182 if (nargout > 3 && iter >= max_iter)
183 exitflag = 0;
184 endif
185 if (nargout > 4)
186 output = struct ("algorithm", "nnls", "iterations", iter);
187 endif
188 if (nargout > 5)
189 lambda = zeros (size (x));
190 lambda(p) = w;
191 endif
192
193endfunction
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)