Mercurial > hg > octave-nkf
annotate scripts/optimization/pqpnonneg.m @ 15144:9cc337ced51a
build: Update OCTAVE_UMFPACK_SEPARATE_SPLIT macro to look for SuiteSparse header file.
* acinclude.m4: Update OCTAVE_UMFPACK_SEPARATE_SPLIT macro to look for
SuiteSparse header file.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 10 Aug 2012 12:56:15 -0700 |
parents | b76f0740940e |
children | bc924baa2c4e |
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 |
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{}) | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10635
diff
changeset
|
28 ## Minimize @code{1/2*x'*c*x + d'*x} subject to @code{@var{x} >= 0}. @var{c} |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10635
diff
changeset
|
29 ## and @var{d} must be real, and @var{c} must be symmetric and positive |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10635
diff
changeset
|
30 ## definite. @var{x0} is an optional initial guess for @var{x}. |
9635 | 31 ## |
32 ## Outputs: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
33 ## |
9635 | 34 ## @itemize @bullet |
35 ## @item minval | |
36 ## | |
37 ## 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
|
38 ## |
9635 | 39 ## @item exitflag |
40 ## | |
41 ## An indicator of convergence. 0 indicates that the iteration count | |
42 ## was exceeded, and therefore convergence was not reached; >0 indicates | |
43 ## that the algorithm converged. (The algorithm is stable and will | |
44 ## converge given enough iterations.) | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
45 ## |
9635 | 46 ## @item output |
47 ## | |
48 ## A structure with two fields: | |
14366
b76f0740940e
doc: Periodic grammar check of documentation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
49 ## |
9635 | 50 ## @itemize @bullet |
51 ## @item "algorithm": The algorithm used ("nnls") | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
52 ## |
9635 | 53 ## @item "iterations": The number of iterations taken. |
54 ## @end itemize | |
10821
693e22af08ae
Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
55 ## |
9635 | 56 ## @item lambda |
57 ## | |
58 ## Not implemented. | |
59 ## @end itemize | |
60 ## @seealso{optimset, lsqnonneg, qp} | |
61 ## @end deftypefn | |
62 | |
13027
b9a89ca0fb75
prevent optimization functions from setting ans in workspace at startup
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
63 ## 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
|
64 ## PKG_ADD: [~] = __all_opts__ ("pqpnonneg"); |
9635 | 65 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
66 ## This is analogical to the lsqnonneg implementation, which is |
9635 | 67 ## implemented from Lawson and Hanson's 1973 algorithm on page |
68 ## 161 of Solving Least Squares Problems. | |
69 ## It shares the convergence guarantees. | |
70 | |
71 function [x, minval, exitflag, output, lambda] = pqpnonneg (c, d, x = [], options = struct ()) | |
72 | |
73 if (nargin == 1 && ischar (c) && strcmp (c, 'defaults')) | |
74 x = optimset ("MaxIter", 1e5); | |
75 return | |
76 endif | |
77 | |
78 if (! (nargin >= 2 && nargin <= 4 && ismatrix (c) && ismatrix (d) && isstruct (options))) | |
79 print_usage (); | |
80 endif | |
81 | |
82 ## Lawson-Hanson Step 1 (LH1): initialize the variables. | |
83 m = rows (c); | |
84 n = columns (c); | |
85 if (m != n) | |
10635
d1978e7364ad
Print name of function in error() string messages.
Rik <octave@nomad.inbox5.com>
parents:
10200
diff
changeset
|
86 error ("pqpnonneg: matrix must be square"); |
9635 | 87 endif |
88 | |
89 if (isempty (x)) | |
90 ## Initial guess is 0s. | |
91 x = zeros (n, 1); | |
92 else | |
93 ## ensure nonnegative guess. | |
94 x = max (x, 0); | |
95 endif | |
96 | |
97 max_iter = optimget (options, "MaxIter", 1e5); | |
98 | |
99 ## Initialize P, according to zero pattern of x. | |
100 p = find (x > 0).'; | |
101 ## Initialize the Cholesky factorization. | |
102 r = chol (c(p, p)); | |
103 usechol = true; | |
104 | |
105 iter = 0; | |
106 | |
107 ## LH3: test for completion. | |
108 while (iter < max_iter) | |
109 while (iter < max_iter) | |
110 iter++; | |
111 | |
112 ## LH6: compute the positive matrix and find the min norm solution | |
113 ## of the positive problem. | |
114 if (usechol) | |
115 xtmp = -(r \ (r' \ d(p))); | |
116 else | |
117 xtmp = -(c(p,p) \ d(p)); | |
118 endif | |
119 idx = find (xtmp < 0); | |
120 | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
121 if (isempty (idx)) |
9635 | 122 ## LH7: tmp solution found, iterate. |
123 x(:) = 0; | |
124 x(p) = xtmp; | |
125 break; | |
126 else | |
127 ## LH8, LH9: find the scaling factor. | |
128 pidx = p(idx); | |
129 sf = x(pidx)./(x(pidx) - xtmp(idx)); | |
130 alpha = min (sf); | |
131 ## LH10: adjust X. | |
132 xx = zeros (n, 1); | |
133 xx(p) = xtmp; | |
134 x += alpha*(xx - x); | |
135 ## LH11: move from P to Z all X == 0. | |
136 ## This corresponds to those indices where minimum of sf is attained. | |
137 idx = idx (sf == alpha); | |
138 p(idx) = []; | |
139 if (usechol) | |
140 ## update the Cholesky factorization. | |
141 r = choldelete (r, idx); | |
142 endif | |
143 endif | |
144 endwhile | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
145 |
9635 | 146 ## compute the gradient. |
147 w = -(d + c*x); | |
148 w(p) = []; | |
149 if (! any (w > 0)) | |
150 if (usechol) | |
151 ## verify the solution achieved using qr updating. | |
152 ## in the best case, this should only take a single step. | |
153 usechol = false; | |
154 continue; | |
155 else | |
156 ## we're finished. | |
157 break; | |
158 endif | |
159 endif | |
160 | |
161 ## find the maximum gradient. | |
162 idx = find (w == max (w)); | |
163 if (numel (idx) > 1) | |
164 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
|
165 "a non-unique solution may be returned due to equal gradients"); |
9635 | 166 idx = idx(1); |
167 endif | |
168 ## move the index from Z to P. Keep P sorted. | |
169 z = [1:n]; z(p) = []; | |
170 zidx = z(idx); | |
171 jdx = 1 + lookup (p, zidx); | |
172 p = [p(1:jdx-1), zidx, p(jdx:end)]; | |
173 if (usechol) | |
174 ## insert the column into the Cholesky factorization. | |
10200
7c1b1c084af1
handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
9635
diff
changeset
|
175 [r, bad] = cholinsert (r, jdx, c(p,zidx)); |
7c1b1c084af1
handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
9635
diff
changeset
|
176 if (bad) |
7c1b1c084af1
handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
9635
diff
changeset
|
177 ## 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
|
178 usechol = false; |
7c1b1c084af1
handle cholesky update failure in pqpnonneg
Jaroslav Hajek <highegg@gmail.com>
parents:
9635
diff
changeset
|
179 endif |
9635 | 180 endif |
181 | |
182 endwhile | |
183 ## LH12: complete. | |
184 | |
185 ## Generate the additional output arguments. | |
186 if (nargout > 1) | |
187 minval = 1/2*(x'*c*x) + d'*x; | |
188 endif | |
189 exitflag = iter; | |
190 if (nargout > 2 && iter >= max_iter) | |
191 exitflag = 0; | |
192 endif | |
193 if (nargout > 3) | |
194 output = struct ("algorithm", "nnls-pqp", "iterations", iter); | |
195 endif | |
196 if (nargout > 4) | |
197 lambda = zeros (size (x)); | |
198 lambda(p) = w; | |
199 endif | |
200 | |
201 endfunction | |
202 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
203 |
9635 | 204 %!test |
205 %! C = [5 2;2 2]; | |
206 %! d = [3; -1]; | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
207 %! assert (pqpnonneg (C, d), [0;0.5], 100*eps); |
9635 | 208 |
209 ## Test equivalence of lsq and pqp | |
210 %!test | |
211 %! C = rand (20, 10); | |
212 %! 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
|
213 %! 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
|
214 |