annotate scripts/optimization/lsqnonneg.m @ 9051:1bf0ce0930be

Grammar check TexInfo in all .m files Cleanup documentation sources to follow a few consistent rules. Spellcheck was NOT done. (but will be in another changeset)
author Rik <rdrider0-list@yahoo.com>
date Fri, 27 Mar 2009 22:31:03 -0700
parents ff61b53eb294
children 36d885c4a1ac
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
3 ## Copyright (C) 2009 VZLU Prague
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
4 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
5 ## This file is part of Octave.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
6 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
7 ## Octave is free software; you can redistribute it and/or modify it
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
8 ## under the terms of the GNU General Public License as published by
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
9 ## the Free Software Foundation; either version 3 of the License, or (at
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
10 ## your option) any later version.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
11 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
12 ## Octave is distributed in the hope that it will be useful, but
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
13 ## WITHOUT ANY WARRANTY; without even the implied warranty of
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
14 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
15 ## General Public License for more details.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
16 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
17 ## You should have received a copy of the GNU General Public License
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
18 ## along with Octave; see the file COPYING. If not, see
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
19 ## <http://www.gnu.org/licenses/>.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
20
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
21 ## -*- texinfo -*-
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
22 ## @deftypefn {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
23 ## @deftypefnx {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}, @var{x0})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
24 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}] =} lsqnonneg (@dots{})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
25 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}] =} lsqnonneg (@dots{})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
26 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}] =} lsqnonneg (@dots{})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
27 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}] =} lsqnonneg (@dots{})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
28 ## @deftypefnx {Function File} {[@var{x}, @var{resnorm}, @var{residual}, @var{exitflag}, @var{output}, @var{lambda}] =} lsqnonneg (@dots{})
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
31 ## initial guess for @var{x}.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
32 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
33 ## Outputs:
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
34 ## @itemize @bullet
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
35 ## @item resnorm
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
36 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
37 ## The squared 2-norm of the residual: norm(@var{c}*@var{x}-@var{d})^2
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
38 ## @item residual
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
39 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
40 ## The residual: @var{d}-@var{c}*@var{x}
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
41 ## @item exitflag
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
42 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
43 ## An indicator of convergence. 0 indicates that the iteration count
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
44 ## was exceeded, and therefore convergence was not reached; >0 indicates
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
45 ## that the algorithm converged. (The algorithm is stable and will
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
46 ## converge given enough iterations.)
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
47 ## @item output
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
48 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
49 ## A structure with two fields:
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
50 ## @itemize @bullet
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
51 ## @item "algorithm": The algorithm used ("nnls")
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
52 ## @item "iterations": The number of iterations taken.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
53 ## @end itemize
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
54 ## @item lambda
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
55 ##
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
56 ## Not implemented.
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
57 ## @end itemize
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
58 ## @seealso{optimset}
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
59 ## @end deftypefn
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
65
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
66 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ())
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
67
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
68 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults'))
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
69 x = optimset ("MaxIter", 1e5);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
70 return
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
71 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
72
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
73 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options)))
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
74 print_usage ();
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
75 endif
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
78 m = rows (c);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
79 n = columns (c);
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
82 x = zeros (n, 1);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
83 else
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
84 ## ensure nonnegative guess.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
85 x = max (x, 0);
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
86 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
87
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
88 useqr = m >= n;
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8407
diff changeset
89 max_iter = optimget (options, "MaxIter", 1e5);
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
90
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
91 ## Initialize P, according to zero pattern of x.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
92 p = find (x > 0).';
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
93 if (useqr)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
94 ## Initialize the QR factorization, economized form.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
98 iter = 0;
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
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
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
101 while (iter < max_iter)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
102 while (iter < max_iter)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
103 iter++;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
104
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
105 ## LH6: compute the positive matrix and find the min norm solution
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
106 ## of the positive problem.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
107 if (useqr)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
108 xtmp = r \ q'*d;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
109 else
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
110 xtmp = c(:,p) \ d;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
111 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
112 idx = find (xtmp < 0);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
113
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
114 if (isempty (idx))
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
115 ## LH7: tmp solution found, iterate.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
116 x(:) = 0;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
117 x(p) = xtmp;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
118 break;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
119 else
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
120 ## LH8, LH9: find the scaling factor.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
121 pidx = p(idx);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
122 sf = x(pidx)./(x(pidx) - xtmp(idx));
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
123 alpha = min (sf);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
124 ## LH10: adjust X.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
125 xx = zeros (n, 1);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
126 xx(p) = xtmp;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
127 x += alpha*(xx - x);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
128 ## LH11: move from P to Z all X == 0.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
129 ## This corresponds to those indices where minimum of sf is attained.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
130 idx = idx (sf == alpha);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
131 p(idx) = [];
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
132 if (useqr)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
133 ## update the QR factorization.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
134 [q, r] = qrdelete (q, r, idx);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
135 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
136 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
137 endwhile
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
138
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
139 ## compute the gradient.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
140 w = c'*(d - c*x);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
141 w(p) = [];
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
142 if (! any (w > 0))
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
143 if (useqr)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
144 ## verify the solution achieved using qr updating.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
145 ## in the best case, this should only take a single step.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
146 useqr = false;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
147 continue;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
148 else
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
149 ## we're finished.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
150 break;
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
151 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
152 endif
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
153
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
154 ## find the maximum gradient.
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
155 idx = find (w == max (w));
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
156 if (numel (idx) > 1)
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
157 warning ("lsqnonneg:nonunique",
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
158 "A non-unique solution may be returned due to equal gradients.");
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
159 idx = idx(1);
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
160 endif
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
161 ## move the index from Z to P. Keep P sorted.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
162 z = [1:n]; z(p) = [];
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
163 zidx = z(idx);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
164 jdx = 1 + lookup (p, zidx);
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
165 p = [p(1:jdx-1), zidx, p(jdx:end)];
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
166 if (useqr)
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
167 ## insert the column into the QR factorization.
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
168 [q, r] = qrinsert (q, r, jdx, c(:,zidx));
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
169 endif
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
170
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
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
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
177 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
178 if (nargout > 2)
8406
0b7566aea627 fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8304
diff changeset
179 residual = d - c*x;
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
180 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
181 exitflag = iter;
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 8407
diff changeset
182 if (nargout > 3 && iter >= max_iter)
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
183 exitflag = 0;
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
184 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
185 if (nargout > 4)
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
186 output = struct ("algorithm", "nnls", "iterations", iter);
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
187 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
188 if (nargout > 5)
8600
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
189 lambda = zeros (size (x));
a6c1aa6f5915 improve lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 8507
diff changeset
190 lambda(p) = w;
7682
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
191 endif
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
192
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
193 endfunction
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
194
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
195 ## Tests
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
196 %!test
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
197 %! C = [1 0;0 1;2 1];
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
198 %! d = [1;3;-2];
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
199 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps)
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
200
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
201 %!test
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
202 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170];
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
203 %! d = [0.8587;0.1781;0.0747;0.8405];
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
204 %! xnew = [0;0.6929];
f7474c83c782 lsqnonneg: new function
bill@denney.ws
parents:
diff changeset
205 %! assert (lsqnonneg (C, d), xnew, 0.0001)