annotate scripts/optimization/pqpnonneg.m @ 20735:418ae0cb752f

Replace ++,-- with in-place operators for performance. * inputParser.m, int2str.m, interp2.m, interp3.m, interpn.m, __isequal__.m, voronoi.m, imshow.m, jet.m, __imread__.m, __imwrite__.m, imwrite_filename.m, dlmwrite.m, importdata.m, strread.m, krylov.m, fact.m, genvarname.m, fminbnd.m, fminsearch.m, fminunc.m, fsolve.m, fzero.m, lsqnonneg.m, pqpnonneg.m, qp.m, sqp.m, annotation.m, legend.m, orient.m, __bar__.m, __contour__.m, __patch__.m, __plt__.m, __scatter__.m, rectangle.m, surfc.m, figure.m, hdl2struct.m, pan.m, __gnuplot_draw_axes__.m, __gnuplot_ginput__.m, rotate3d.m, struct2hdl.m, zoom.m, freqz.m, periodogram.m, gmres.m, pcg.m, pcr.m, __sprand__.m, sprandsym.m, treelayout.m, treeplot.m, hadamard.m, hilb.m, iqr.m, assert.m, test.m: Replace ++,-- with in-place operators for performance.
author Rik <rik@octave.org>
date Sun, 22 Nov 2015 20:25:17 -0800
parents f1d0f506ee78
children 516bb87ea72e
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 17744
diff changeset
1 ## Copyright (C) 2008-2015 Bill Denney
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
2 ## Copyright (C) 2008 Jaroslav Hajek
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
3 ## Copyright (C) 2009 VZLU Prague
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
4 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
5 ## This file is part of Octave.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
6 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
7 ## Octave is free software; you can redistribute it and/or modify it
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
8 ## under the terms of the GNU General Public License as published by
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
9 ## the Free Software Foundation; either version 3 of the License, or (at
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
10 ## your option) any later version.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
11 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
12 ## Octave is distributed in the hope that it will be useful, but
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
13 ## WITHOUT ANY WARRANTY; without even the implied warranty of
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
14 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
15 ## General Public License for more details.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
16 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
17 ## You should have received a copy of the GNU General Public License
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
18 ## along with Octave; see the file COPYING. If not, see
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
19 ## <http://www.gnu.org/licenses/>.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
20
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
21 ## -*- texinfo -*-
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10635
diff changeset
22 ## @deftypefn {Function File} {@var{x} =} pqpnonneg (@var{c}, @var{d})
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
23 ## @deftypefnx {Function File} {@var{x} =} pqpnonneg (@var{c}, @var{d}, @var{x0})
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
24 ## @deftypefnx {Function File} {[@var{x}, @var{minval}] =} pqpnonneg (@dots{})
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
25 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}] =} pqpnonneg (@dots{})
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
26 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}, @var{output}] =} pqpnonneg (@dots{})
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
27 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}, @var{output}, @var{lambda}] =} pqpnonneg (@dots{})
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
28 ## Minimize @code{1/2*x'*c*x + d'*x} subject to @code{@var{x} >= 0}.
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
29 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
30 ## @var{c} ## and @var{d} must be real, and @var{c} must be symmetric and
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
31 ## positive definite.
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
32 ##
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
33 ## @var{x0} is an optional initial guess for @var{x}.
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
34 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
35 ## Outputs:
14366
b76f0740940e doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
36 ##
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
37 ## @itemize @bullet
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
38 ## @item minval
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
39 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
40 ## The minimum attained model value, 1/2*xmin'*c*xmin + d'*xmin
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
41 ##
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
42 ## @item exitflag
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
43 ##
20165
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
44 ## An indicator of convergence. 0 indicates that the iteration count was
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
45 ## exceeded, and therefore convergence was not reached; >0 indicates that the
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
46 ## algorithm converged. (The algorithm is stable and will converge given
f1d0f506ee78 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
47 ## enough iterations.)
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
48 ##
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
49 ## @item output
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
50 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
51 ## A structure with two fields:
14366
b76f0740940e doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
52 ##
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
53 ## @itemize @bullet
17281
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 14366
diff changeset
54 ## @item @qcode{"algorithm"}: The algorithm used (@qcode{"nnls"})
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
55 ##
17281
bc924baa2c4e doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents: 14366
diff changeset
56 ## @item @qcode{"iterations"}: The number of iterations taken.
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
57 ## @end itemize
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10793
diff changeset
58 ##
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
59 ## @item lambda
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
60 ##
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
61 ## Not implemented.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
62 ## @end itemize
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
63 ## @seealso{optimset, lsqnonneg, qp}
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
64 ## @end deftypefn
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
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__ ("pqpnonneg");
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
68
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
69 ## This is analogical to the lsqnonneg implementation, which is
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
70 ## implemented from Lawson and Hanson's 1973 algorithm on page
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
71 ## 161 of Solving Least Squares Problems.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
72 ## It shares the convergence guarantees.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
73
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
74 function [x, minval, exitflag, output, lambda] = pqpnonneg (c, d, x = [], options = struct ())
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
75
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
76 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults'))
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
77 x = optimset ("MaxIter", 1e5);
17312
088d014a7fe2 Use semicolon after "return" statement in core m-files.
Rik <rik@octave.org>
parents: 17281
diff changeset
78 return;
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
79 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
80
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
81 if (nargin < 2 || nargin > 4
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
82 || ! (ismatrix (c) && ismatrix (d) && isstruct (options)))
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
83 print_usage ();
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
84 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
85
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
86 ## Lawson-Hanson Step 1 (LH1): initialize the variables.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
87 m = rows (c);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
88 n = columns (c);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
89 if (m != n)
10635
d1978e7364ad Print name of function in error() string messages.
Rik <octave@nomad.inbox5.com>
parents: 10200
diff changeset
90 error ("pqpnonneg: matrix must be square");
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
91 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
92
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
93 if (isempty (x))
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
94 ## Initial guess is 0s.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
95 x = zeros (n, 1);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
96 else
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
97 ## ensure nonnegative guess.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
98 x = max (x, 0);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
99 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
100
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
101 max_iter = optimget (options, "MaxIter", 1e5);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
102
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
103 ## Initialize P, according to zero pattern of x.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
104 p = find (x > 0).';
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
105 ## Initialize the Cholesky factorization.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
106 r = chol (c(p, p));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
107 usechol = true;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
108
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
109 iter = 0;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
110
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
111 ## LH3: test for completion.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
112 while (iter < max_iter)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
113 while (iter < max_iter)
20735
418ae0cb752f Replace ++,-- with in-place operators for performance.
Rik <rik@octave.org>
parents: 20165
diff changeset
114 iter += 1;
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
115
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
116 ## LH6: compute the positive matrix and find the min norm solution
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
117 ## of the positive problem.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
118 if (usechol)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
119 xtmp = -(r \ (r' \ d(p)));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
120 else
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
121 xtmp = -(c(p,p) \ d(p));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
122 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
123 idx = find (xtmp < 0);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
124
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
125 if (isempty (idx))
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
126 ## LH7: tmp solution found, iterate.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
127 x(:) = 0;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
128 x(p) = xtmp;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
129 break;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
130 else
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
131 ## LH8, LH9: find the scaling factor.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
132 pidx = p(idx);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
133 sf = x(pidx)./(x(pidx) - xtmp(idx));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
134 alpha = min (sf);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
135 ## LH10: adjust X.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
136 xx = zeros (n, 1);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
137 xx(p) = xtmp;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
138 x += alpha*(xx - x);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
139 ## LH11: move from P to Z all X == 0.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
140 ## This corresponds to those indices where minimum of sf is attained.
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
141 idx = idx(sf == alpha);
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
142 p(idx) = [];
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
143 if (usechol)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
144 ## update the Cholesky factorization.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
145 r = choldelete (r, idx);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
146 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
147 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
148 endwhile
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
149
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
150 ## compute the gradient.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
151 w = -(d + c*x);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
152 w(p) = [];
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
153 if (! any (w > 0))
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
154 if (usechol)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
155 ## verify the solution achieved using qr updating.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
156 ## in the best case, this should only take a single step.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
157 usechol = false;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
158 continue;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
159 else
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
160 ## we're finished.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
161 break;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
162 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
163 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
164
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
165 ## find the maximum gradient.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
166 idx = find (w == max (w));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
167 if (numel (idx) > 1)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
168 warning ("pqpnonneg:nonunique",
11588
d5bd2766c640 style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
169 "a non-unique solution may be returned due to equal gradients");
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
170 idx = idx(1);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
171 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
172 ## move the index from Z to P. Keep P sorted.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
173 z = [1:n]; z(p) = [];
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
174 zidx = z(idx);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
175 jdx = 1 + lookup (p, zidx);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
176 p = [p(1:jdx-1), zidx, p(jdx:end)];
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
177 if (usechol)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
178 ## insert the column into the Cholesky factorization.
10200
7c1b1c084af1 handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 9635
diff changeset
179 [r, bad] = cholinsert (r, jdx, c(p,zidx));
7c1b1c084af1 handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 9635
diff changeset
180 if (bad)
7c1b1c084af1 handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 9635
diff changeset
181 ## If the insertion failed, we switch off updates and go on.
7c1b1c084af1 handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 9635
diff changeset
182 usechol = false;
7c1b1c084af1 handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents: 9635
diff changeset
183 endif
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
184 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
185
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
186 endwhile
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
187 ## LH12: complete.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
188
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
189 ## Generate the additional output arguments.
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
190 if (nargout > 1)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
191 minval = 1/2*(x'*c*x) + d'*x;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
192 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
193 exitflag = iter;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
194 if (nargout > 2 && iter >= max_iter)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
195 exitflag = 0;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
196 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
197 if (nargout > 3)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
198 output = struct ("algorithm", "nnls-pqp", "iterations", iter);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
199 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
200 if (nargout > 4)
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
201 lambda = zeros (size (x));
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
202 lambda(p) = w;
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
203 endif
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
204
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
205 endfunction
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
206
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
207
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
208 %!test
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
209 %! C = [5 2;2 2];
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
210 %! d = [3; -1];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
211 %! assert (pqpnonneg (C, d), [0;0.5], 100*eps);
9635
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
212
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
213 ## Test equivalence of lsq and pqp
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
214 %!test
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
215 %! C = rand (20, 10);
36d885c4a1ac implement pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
diff changeset
216 %! d = rand (20, 1);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
217 %! assert (pqpnonneg (C'*C, -C'*d), lsqnonneg (C, d), 100*eps);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
218