1## Copyright (C) 2008 Bill Denney
2## Copyright (C) 2008 Jaroslav Hajek
3## Copyright (C) 2009 VZLU Prague
5## This file is part of Octave.
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.
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.
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/>.
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}.
37## The squared 2-norm of the residual: norm(@var{c}*@var{x}-@var{d})^2
40## The residual: @var{d}-@var{c}*@var{x}
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.)
49## A structure with two fields:
51## @item "algorithm": The algorithm used ("nnls")
52## @item "iterations": The number of iterations taken.
58## @seealso{optimset, pqpnonneg}
61## PKG_ADD: __all_opts__ ("lsqnonneg");
63## This is implemented from Lawson and Hanson's 1973 algorithm on page
64## 161 of Solving Least Squares Problems.
66function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ())
68 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults'))
69 x = optimset ("MaxIter", 1e5);
73 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options)))
77 ## Lawson-Hanson Step 1 (LH1): initialize the variables.
81 ## Initial guess is 0s.
84 ## ensure nonnegative guess.
89 max_iter = optimget (options, "MaxIter", 1e5);
91 ## Initialize P, according to zero pattern of x.
94 ## Initialize the QR factorization, economized form.
95 [q, r] = qr (c(:,p), 0);
100 ## LH3: test for completion.
101 while (iter < max_iter)
102 while (iter < max_iter)
105 ## LH6: compute the positive matrix and find the min norm solution
106 ## of the positive problem.
112 idx = find (xtmp < 0);
115 ## LH7: tmp solution found, iterate.
120 ## LH8, LH9: find the scaling factor.
122 sf = x(pidx)./(x(pidx) - xtmp(idx));
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);
133 ## update the QR factorization.
134 [q, r] = qrdelete (q, r, idx);
139 ## compute the gradient.
144 ## verify the solution achieved using qr updating.
145 ## in the best case, this should only take a single step.
154 ## find the maximum gradient.
155 idx = find (w == max (w));
157 warning ("lsqnonneg:nonunique",
158 "A non-unique solution may be returned due to equal gradients.");
161 ## move the index from Z to P. Keep P sorted.
162 z = [1:n]; z(p) = [];
164 jdx = 1 + lookup (p, zidx);
165 p = [p(1:jdx-1), zidx, p(jdx:end)];
167 ## insert the column into the QR factorization.
168 [q, r] = qrinsert (q, r, jdx, c(:,zidx));
174 ## Generate the additional output arguments.
176 resnorm = norm (c*x - d) ^ 2;
182 if (nargout > 3 && iter >= max_iter)
186 output = struct ("algorithm", "nnls", "iterations", iter);
189 lambda = zeros (size (x));
199%! assert (lsqnonneg (C, d), [0;0.5], 100*eps)
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];
205%! assert (lsqnonneg (C, d), xnew, 0.0001)