Mercurial > octave-nkf
annotate scripts/optimization/lsqnonneg.m @ 14363:f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
* wavread.m, acosd.m, acot.m, acotd.m, acoth.m, acsc.m, acscd.m, acsch.m,
asec.m, asecd.m, asech.m, asind.m, atand.m, cosd.m, cot.m, cotd.m, coth.m,
csc.m, cscd.m, csch.m, sec.m, secd.m, sech.m, sind.m, tand.m, accumarray.m,
accumdim.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, cart2pol.m, cart2sph.m,
celldisp.m, chop.m, circshift.m, colon.m, common_size.m, cplxpair.m,
cumtrapz.m, curl.m, dblquad.m, deal.m, divergence.m, flipdim.m, fliplr.m,
flipud.m, genvarname.m, gradient.m, idivide.m, int2str.m, interp1.m,
interp1q.m, interp2.m, interp3.m, interpft.m, interpn.m, isa.m, isdir.m,
isequal.m, isequalwithequalnans.m, issquare.m, logspace.m, nargchk.m,
narginchk.m, nargoutchk.m, nextpow2.m, nthargout.m, num2str.m, pol2cart.m,
polyarea.m, postpad.m, prepad.m, profile.m, profshow.m, quadgk.m, quadv.m,
randi.m, rat.m, repmat.m, rot90.m, rotdim.m, shift.m, shiftdim.m, sph2cart.m,
structfun.m, trapz.m, triplequad.m, convhull.m, dsearch.m, dsearchn.m,
griddata3.m, griddatan.m, rectint.m, tsearchn.m, __makeinfo__.m, doc.m,
get_first_help_sentence.m, help.m, type.m, unimplemented.m, which.m, imread.m,
imwrite.m, dlmwrite.m, fileread.m, is_valid_file_id.m, strread.m, textread.m,
textscan.m, commutation_matrix.m, cond.m, condest.m, cross.m,
duplication_matrix.m, expm.m, housh.m, isdefinite.m, ishermitian.m,
issymmetric.m, logm.m, normest.m, null.m, onenormest.m, orth.m, planerot.m,
qzhess.m, rank.m, rref.m, trace.m, vech.m, ans.m, bincoeff.m, bug_report.m,
bzip2.m, comma.m, compare_versions.m, computer.m, edit.m, fileparts.m,
fullfile.m, getfield.m, gzip.m, info.m, inputname.m, isappdata.m, isdeployed.m,
ismac.m, ispc.m, isunix.m, list_primes.m, ls.m, mexext.m, namelengthmax.m,
news.m, orderfields.m, paren.m, recycle.m, rmappdata.m, semicolon.m,
setappdata.m, setfield.m, substruct.m, symvar.m, ver.m, version.m,
warning_ids.m, xor.m, fminbnd.m, fsolve.m, fzero.m, lsqnonneg.m, optimset.m,
pqpnonneg.m, sqp.m, matlabroot.m, __gnuplot_drawnow__.m,
__plt_get_axis_arg__.m, ancestor.m, cla.m, clf.m, close.m, colorbar.m,
colstyle.m, comet3.m, contourc.m, figure.m, gca.m, gcbf.m, gcbo.m, gcf.m,
ginput.m, graphics_toolkit.m, gtext.m, hggroup.m, hist.m, hold.m, isfigure.m,
ishghandle.m, ishold.m, isocolors.m, isonormals.m, isosurface.m, isprop.m,
legend.m, line.m, loglog.m, loglogerr.m, meshgrid.m, ndgrid.m, newplot.m,
orient.m, patch.m, plot3.m, plotyy.m, __print_parse_opts__.m, quiver3.m,
refreshdata.m, ribbon.m, semilogx.m, semilogxerr.m, semilogy.m, stem.m,
stem3.m, subplot.m, title.m, uigetfile.m, view.m, whitebg.m, compan.m, conv.m,
deconv.m, mkpp.m, mpoles.m, pchip.m, poly.m, polyaffine.m, polyder.m,
polyfit.m, polygcd.m, polyint.m, polyout.m, polyval.m, polyvalm.m, ppder.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, intersect.m,
ismember.m, powerset.m, setdiff.m, setxor.m, union.m, unique.m,
autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, fftconv.m, fftfilt.m,
fftshift.m, freqz.m, hamming.m, hanning.m, ifftshift.m, sinc.m, sinetone.m,
sinewave.m, unwrap.m, bicg.m, bicgstab.m, gmres.m, gplot.m, nonzeros.m, pcg.m,
pcr.m, spaugment.m, spconvert.m, spdiags.m, speye.m, spfun.m, spones.m,
sprand.m, sprandsym.m, spstats.m, spy.m, svds.m, treelayout.m, bessel.m,
beta.m, betaln.m, factor.m, factorial.m, isprime.m, lcm.m, legendre.m,
nchoosek.m, nthroot.m, perms.m, pow2.m, primes.m, reallog.m, realpow.m,
realsqrt.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, rosser.m,
vander.m, __finish__.m, center.m, cloglog.m, corr.m, cov.m, gls.m, histc.m,
iqr.m, kendall.m, kurtosis.m, logit.m, mahalanobis.m, mean.m, meansq.m,
median.m, mode.m, moment.m, ols.m, ppplot.m, prctile.m, probit.m, quantile.m,
range.m, ranks.m, run_count.m, runlength.m, skewness.m, spearman.m,
statistics.m, std.m, table.m, var.m, zscore.m, betacdf.m, betainv.m, betapdf.m,
betarnd.m, binocdf.m, binoinv.m, binopdf.m, binornd.m, cauchy_cdf.m,
cauchy_inv.m, cauchy_pdf.m, cauchy_rnd.m, chi2cdf.m, chi2inv.m, chi2pdf.m,
chi2rnd.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m,
empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, expcdf.m,
expinv.m, exppdf.m, exprnd.m, fcdf.m, finv.m, fpdf.m, frnd.m, gamcdf.m,
gaminv.m, gampdf.m, gamrnd.m, geocdf.m, geoinv.m, geopdf.m, geornd.m,
hygecdf.m, hygeinv.m, hygepdf.m, hygernd.m, kolmogorov_smirnov_cdf.m,
laplace_cdf.m, laplace_inv.m, laplace_pdf.m, laplace_rnd.m, logistic_cdf.m,
logistic_inv.m, logistic_pdf.m, logistic_rnd.m, logncdf.m, logninv.m,
lognpdf.m, lognrnd.m, nbincdf.m, nbininv.m, nbinpdf.m, nbinrnd.m, normcdf.m,
norminv.m, normpdf.m, normrnd.m, poisscdf.m, poissinv.m, poisspdf.m,
poissrnd.m, stdnormal_cdf.m, stdnormal_inv.m, stdnormal_pdf.m, stdnormal_rnd.m,
tcdf.m, tinv.m, tpdf.m, trnd.m, unidcdf.m, unidinv.m, unidpdf.m, unidrnd.m,
unifcdf.m, unifinv.m, unifpdf.m, unifrnd.m, wblcdf.m, wblinv.m, wblpdf.m,
wblrnd.m, kolmogorov_smirnov_test.m, kruskal_wallis_test.m, base2dec.m,
bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m,
findstr.m, hex2dec.m, index.m, isletter.m, mat2str.m, rindex.m, str2num.m,
strcat.m, strjust.m, strmatch.m, strsplit.m, strtok.m, strtrim.m, strtrunc.m,
substr.m, validatestring.m, demo.m, example.m, fail.m, speed.m, addtodate.m,
asctime.m, clock.m, ctime.m, date.m, datenum.m, datetick.m, datevec.m,
eomday.m, etime.m, is_leap_year.m, now.m:
Use Octave coding conventions in all m-file %!test blocks
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Mon, 13 Feb 2012 07:29:44 -0800 |
parents | 72c96de7a403 |
children | b76f0740940e |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13027
diff
changeset
|
1 ## Copyright (C) 2008-2012 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 -*- | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
9635
diff
changeset
|
22 ## @deftypefn {Function File} {@var{x} =} lsqnonneg (@var{c}, @var{d}) |
7682 | 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{}) | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
29 ## Minimize @code{norm (@var{c}*@var{x} - d)} subject to |
10846
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
30 ## @code{@var{x} >= 0}. @var{c} and @var{d} must be real. @var{x0} is an |
a4f482e66b65
Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents:
10821
diff
changeset
|
31 ## optional initial guess for @var{x}. |
7682 | 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 | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
38 ## |
7682 | 39 ## @item residual |
40 ## | |
41 ## The residual: @var{d}-@var{c}*@var{x} | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
42 ## |
7682 | 43 ## @item exitflag |
44 ## | |
45 ## An indicator of convergence. 0 indicates that the iteration count | |
46 ## was exceeded, and therefore convergence was not reached; >0 indicates | |
47 ## that the algorithm converged. (The algorithm is stable and will | |
48 ## converge given enough iterations.) | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
49 ## |
7682 | 50 ## @item output |
51 ## | |
52 ## A structure with two fields: | |
53 ## @itemize @bullet | |
54 ## @item "algorithm": The algorithm used ("nnls") | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## |
7682 | 56 ## @item "iterations": The number of iterations taken. |
57 ## @end itemize | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
58 ## |
7682 | 59 ## @item lambda |
60 ## | |
61 ## Not implemented. | |
62 ## @end itemize | |
9635 | 63 ## @seealso{optimset, pqpnonneg} |
7682 | 64 ## @end deftypefn |
65 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
66 ## PKG_ADD: ## Discard result to avoid polluting workspace with ans at startup. |
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
67 ## PKG_ADD: [~] = __all_opts__ ("lsqnonneg"); |
8648
ff61b53eb294
optimization: use PKG_ADD: comments instead of PKG_ADD file
John W. Eaton <jwe@octave.org>
parents:
8600
diff
changeset
|
68 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
69 ## 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
|
70 ## 161 of Solving Least Squares Problems. |
7682 | 71 |
8600 | 72 function [x, resnorm, residual, exitflag, output, lambda] = lsqnonneg (c, d, x = [], options = struct ()) |
73 | |
74 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
75 x = optimset ("MaxIter", 1e5); | |
76 return | |
77 endif | |
78 | |
79 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options))) | |
80 print_usage (); | |
81 endif | |
7682 | 82 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
83 ## Lawson-Hanson Step 1 (LH1): initialize the variables. |
8600 | 84 m = rows (c); |
85 n = columns (c); | |
7682 | 86 if (isempty (x)) |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
87 ## Initial guess is 0s. |
8600 | 88 x = zeros (n, 1); |
89 else | |
90 ## ensure nonnegative guess. | |
91 x = max (x, 0); | |
7682 | 92 endif |
93 | |
8600 | 94 useqr = m >= n; |
8507 | 95 max_iter = optimget (options, "MaxIter", 1e5); |
7682 | 96 |
8600 | 97 ## Initialize P, according to zero pattern of x. |
98 p = find (x > 0).'; | |
99 if (useqr) | |
100 ## Initialize the QR factorization, economized form. | |
101 [q, r] = qr (c(:,p), 0); | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
102 endif |
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
103 |
7682 | 104 iter = 0; |
8600 | 105 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
106 ## LH3: test for completion. |
8600 | 107 while (iter < max_iter) |
108 while (iter < max_iter) | |
109 iter++; | |
110 | |
111 ## LH6: compute the positive matrix and find the min norm solution | |
112 ## of the positive problem. | |
113 if (useqr) | |
114 xtmp = r \ q'*d; | |
115 else | |
116 xtmp = c(:,p) \ d; | |
117 endif | |
118 idx = find (xtmp < 0); | |
119 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
120 if (isempty (idx)) |
8600 | 121 ## LH7: tmp solution found, iterate. |
122 x(:) = 0; | |
123 x(p) = xtmp; | |
124 break; | |
125 else | |
126 ## LH8, LH9: find the scaling factor. | |
127 pidx = p(idx); | |
128 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
129 alpha = min (sf); | |
130 ## LH10: adjust X. | |
131 xx = zeros (n, 1); | |
132 xx(p) = xtmp; | |
133 x += alpha*(xx - x); | |
134 ## LH11: move from P to Z all X == 0. | |
135 ## This corresponds to those indices where minimum of sf is attained. | |
136 idx = idx (sf == alpha); | |
137 p(idx) = []; | |
138 if (useqr) | |
139 ## update the QR factorization. | |
140 [q, r] = qrdelete (q, r, idx); | |
141 endif | |
142 endif | |
143 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
144 |
8600 | 145 ## compute the gradient. |
146 w = c'*(d - c*x); | |
147 w(p) = []; | |
148 if (! any (w > 0)) | |
149 if (useqr) | |
150 ## verify the solution achieved using qr updating. | |
151 ## in the best case, this should only take a single step. | |
152 useqr = false; | |
153 continue; | |
154 else | |
155 ## we're finished. | |
156 break; | |
157 endif | |
158 endif | |
159 | |
160 ## find the maximum gradient. | |
7682 | 161 idx = find (w == max (w)); |
162 if (numel (idx) > 1) | |
163 warning ("lsqnonneg:nonunique", | |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
164 "a non-unique solution may be returned due to equal gradients"); |
7682 | 165 idx = idx(1); |
166 endif | |
8600 | 167 ## move the index from Z to P. Keep P sorted. |
168 z = [1:n]; z(p) = []; | |
169 zidx = z(idx); | |
170 jdx = 1 + lookup (p, zidx); | |
171 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
172 if (useqr) | |
173 ## insert the column into the QR factorization. | |
174 [q, r] = qrinsert (q, r, jdx, c(:,zidx)); | |
175 endif | |
7682 | 176 |
177 endwhile | |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
178 ## LH12: complete. |
7682 | 179 |
7697
0bdfff62cc49
lsqnonneg: use optimset, correctly index Z and P in main loop
bill@denney.ws
parents:
7682
diff
changeset
|
180 ## Generate the additional output arguments. |
7682 | 181 if (nargout > 1) |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
182 resnorm = norm (c*x - d) ^ 2; |
7682 | 183 endif |
184 if (nargout > 2) | |
8406
0b7566aea627
fix & optimize lsqnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
8304
diff
changeset
|
185 residual = d - c*x; |
7682 | 186 endif |
187 exitflag = iter; | |
8507 | 188 if (nargout > 3 && iter >= max_iter) |
7682 | 189 exitflag = 0; |
190 endif | |
191 if (nargout > 4) | |
192 output = struct ("algorithm", "nnls", "iterations", iter); | |
193 endif | |
194 if (nargout > 5) | |
8600 | 195 lambda = zeros (size (x)); |
196 lambda(p) = w; | |
7682 | 197 endif |
198 | |
199 endfunction | |
200 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
201 |
7682 | 202 %!test |
203 %! C = [1 0;0 1;2 1]; | |
204 %! d = [1;3;-2]; | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
205 %! assert (lsqnonneg (C, d), [0;0.5], 100*eps); |
7682 | 206 |
207 %!test | |
208 %! C = [0.0372 0.2869;0.6861 0.7071;0.6233 0.6245;0.6344 0.6170]; | |
209 %! d = [0.8587;0.1781;0.0747;0.8405]; | |
210 %! xnew = [0;0.6929]; | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
211 %! assert (lsqnonneg (C, d), xnew, 0.0001); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
212 |