Mercurial > octave-nkf
annotate scripts/optimization/pqpnonneg.m @ 20654:b65888ec820e draft default tip gccjit
dmalcom gcc jit import
author | Stefan Mahr <dac922@gmx.de> |
---|---|
date | Fri, 27 Feb 2015 16:59:36 +0100 |
parents | f1d0f506ee78 |
children |
rev | line source |
---|---|
19731
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 | 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 -*- | |
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 | 23 ## @deftypefnx {Function File} {@var{x} =} pqpnonneg (@var{c}, @var{d}, @var{x0}) |
24 ## @deftypefnx {Function File} {[@var{x}, @var{minval}] =} pqpnonneg (@dots{}) | |
25 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}] =} pqpnonneg (@dots{}) | |
26 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}, @var{output}] =} pqpnonneg (@dots{}) | |
27 ## @deftypefnx {Function File} {[@var{x}, @var{minval}, @var{exitflag}, @var{output}, @var{lambda}] =} pqpnonneg (@dots{}) | |
20200
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19867
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:
19867
diff
changeset
|
29 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19867
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:
19867
diff
changeset
|
31 ## positive definite. |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19867
diff
changeset
|
32 ## |
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19867
diff
changeset
|
33 ## @var{x0} is an optional initial guess for @var{x}. |
9635 | 34 ## |
35 ## Outputs: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
36 ## |
9635 | 37 ## @itemize @bullet |
38 ## @item minval | |
39 ## | |
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 | 42 ## @item exitflag |
43 ## | |
20200
f1d0f506ee78
doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19867
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:
19867
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:
19867
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:
19867
diff
changeset
|
47 ## enough iterations.) |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
48 ## |
9635 | 49 ## @item output |
50 ## | |
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 | 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 | 57 ## @end itemize |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
58 ## |
9635 | 59 ## @item lambda |
60 ## | |
61 ## Not implemented. | |
62 ## @end itemize | |
63 ## @seealso{optimset, lsqnonneg, qp} | |
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__ ("pqpnonneg"); |
9635 | 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 | 70 ## implemented from Lawson and Hanson's 1973 algorithm on page |
71 ## 161 of Solving Least Squares Problems. | |
72 ## It shares the convergence guarantees. | |
73 | |
74 function [x, minval, exitflag, output, lambda] = pqpnonneg (c, d, x = [], options = struct ()) | |
75 | |
76 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
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 | 79 endif |
80 | |
19867
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
81 if (nargin < 2 || nargin > 4 |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
82 || ! (ismatrix (c) && ismatrix (d) && isstruct (options))) |
9635 | 83 print_usage (); |
84 endif | |
85 | |
86 ## Lawson-Hanson Step 1 (LH1): initialize the variables. | |
87 m = rows (c); | |
88 n = columns (c); | |
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 | 91 endif |
92 | |
93 if (isempty (x)) | |
94 ## Initial guess is 0s. | |
95 x = zeros (n, 1); | |
96 else | |
97 ## ensure nonnegative guess. | |
98 x = max (x, 0); | |
99 endif | |
100 | |
101 max_iter = optimget (options, "MaxIter", 1e5); | |
102 | |
103 ## Initialize P, according to zero pattern of x. | |
104 p = find (x > 0).'; | |
105 ## Initialize the Cholesky factorization. | |
106 r = chol (c(p, p)); | |
107 usechol = true; | |
108 | |
109 iter = 0; | |
110 | |
111 ## LH3: test for completion. | |
112 while (iter < max_iter) | |
113 while (iter < max_iter) | |
114 iter++; | |
115 | |
116 ## LH6: compute the positive matrix and find the min norm solution | |
117 ## of the positive problem. | |
118 if (usechol) | |
119 xtmp = -(r \ (r' \ d(p))); | |
120 else | |
121 xtmp = -(c(p,p) \ d(p)); | |
122 endif | |
123 idx = find (xtmp < 0); | |
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 | 126 ## LH7: tmp solution found, iterate. |
127 x(:) = 0; | |
128 x(p) = xtmp; | |
129 break; | |
130 else | |
131 ## LH8, LH9: find the scaling factor. | |
132 pidx = p(idx); | |
133 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
134 alpha = min (sf); | |
135 ## LH10: adjust X. | |
136 xx = zeros (n, 1); | |
137 xx(p) = xtmp; | |
138 x += alpha*(xx - x); | |
139 ## LH11: move from P to Z all X == 0. | |
140 ## This corresponds to those indices where minimum of sf is attained. | |
19867
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19731
diff
changeset
|
141 idx = idx(sf == alpha); |
9635 | 142 p(idx) = []; |
143 if (usechol) | |
144 ## update the Cholesky factorization. | |
145 r = choldelete (r, idx); | |
146 endif | |
147 endif | |
148 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
149 |
9635 | 150 ## compute the gradient. |
151 w = -(d + c*x); | |
152 w(p) = []; | |
153 if (! any (w > 0)) | |
154 if (usechol) | |
155 ## verify the solution achieved using qr updating. | |
156 ## in the best case, this should only take a single step. | |
157 usechol = false; | |
158 continue; | |
159 else | |
160 ## we're finished. | |
161 break; | |
162 endif | |
163 endif | |
164 | |
165 ## find the maximum gradient. | |
166 idx = find (w == max (w)); | |
167 if (numel (idx) > 1) | |
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 | 170 idx = idx(1); |
171 endif | |
172 ## move the index from Z to P. Keep P sorted. | |
173 z = [1:n]; z(p) = []; | |
174 zidx = z(idx); | |
175 jdx = 1 + lookup (p, zidx); | |
176 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
177 if (usechol) | |
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 | 184 endif |
185 | |
186 endwhile | |
187 ## LH12: complete. | |
188 | |
189 ## Generate the additional output arguments. | |
190 if (nargout > 1) | |
191 minval = 1/2*(x'*c*x) + d'*x; | |
192 endif | |
193 exitflag = iter; | |
194 if (nargout > 2 && iter >= max_iter) | |
195 exitflag = 0; | |
196 endif | |
197 if (nargout > 3) | |
198 output = struct ("algorithm", "nnls-pqp", "iterations", iter); | |
199 endif | |
200 if (nargout > 4) | |
201 lambda = zeros (size (x)); | |
202 lambda(p) = w; | |
203 endif | |
204 | |
205 endfunction | |
206 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
207 |
9635 | 208 %!test |
209 %! C = [5 2;2 2]; | |
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 | 212 |
213 ## Test equivalence of lsq and pqp | |
214 %!test | |
215 %! C = rand (20, 10); | |
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 |