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