1 ## Copyright (C) 2006-2012 David Bateman
3 ## This file is part of Octave.
5 ## Octave is free software; you can redistribute it and/or modify it
6 ## under the terms of the GNU General Public License as published by
7 ## the Free Software Foundation; either version 3 of the License, or (at
8 ## your option) any later version.
10 ## Octave is distributed in the hope that it will be useful, but
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 ## General Public License for more details.
15 ## You should have received a copy of the GNU General Public License
16 ## along with Octave; see the file COPYING. If not, see
17 ## <http://www.gnu.org/licenses/>.
20 ## @deftypefn {Function File} {@var{s} =} svds (@var{A})
21 ## @deftypefnx {Function File} {@var{s} =} svds (@var{A}, @var{k})
22 ## @deftypefnx {Function File} {@var{s} =} svds (@var{A}, @var{k}, @var{sigma})
23 ## @deftypefnx {Function File} {@var{s} =} svds (@var{A}, @var{k}, @var{sigma}, @var{opts})
24 ## @deftypefnx {Function File} {[@var{u}, @var{s}, @var{v}] =} svds (@dots{})
25 ## @deftypefnx {Function File} {[@var{u}, @var{s}, @var{v}, @var{flag}] =} svds (@dots{})
27 ## Find a few singular values of the matrix @var{A}. The singular values
28 ## are calculated using
32 ## [@var{m}, @var{n}] = size (@var{A});
33 ## @var{s} = eigs ([sparse(@var{m}, @var{m}), @var{A};
34 ## @var{A}', sparse(@var{n}, @var{n})])
38 ## The eigenvalues returned by @code{eigs} correspond to the singular values
39 ## of @var{A}. The number of singular values to calculate is given by @var{k}
42 ## The argument @var{sigma} specifies which singular values to find. When
43 ## @var{sigma} is the string 'L', the default, the largest singular values of
44 ## @var{A} are found. Otherwise, @var{sigma} must be a real scalar and the
45 ## singular values closest to @var{sigma} are found. As a corollary,
46 ## @code{@var{sigma} = 0} finds the smallest singular values. Note that for
47 ## relatively small values of @var{sigma}, there is a chance that the requested
48 ## number of singular values will not be found. In that case @var{sigma}
49 ## should be increased.
51 ## @var{opts} is a structure defining options that @code{svds} will pass
52 ## to @code{eigs}. The possible fields of this structure are documented in
53 ## @code{eigs}. By default, @code{svds} sets the following three fields:
57 ## The required convergence tolerance for the singular values. The default
58 ## value is 1e-10. @code{eigs} is passed @code{@var{tol} / sqrt(2)}.
61 ## The maximum number of iterations. The default is 300.
64 ## The level of diagnostic printout (0|1|2). If @code{disp} is 0 then
65 ## diagnostics are disabled. The default value is 0.
68 ## If more than one output is requested then @code{svds} will return an
69 ## approximation of the singular value decomposition of @var{A}
72 ## @var{A}_approx = @var{u}*@var{s}*@var{v}'
76 ## where @var{A}_approx is a matrix of size @var{A} but only rank @var{k}.
78 ## @var{flag} returns 0 if the algorithm has succesfully converged, and 1
79 ## otherwise. The test for convergence is
83 ## norm (@var{A}*@var{v} - @var{u}*@var{s}, 1) <= @var{tol} * norm (@var{A}, 1)
87 ## @code{svds} is best for finding only a few singular values from a large
88 ## sparse matrix. Otherwise, @code{svd (full(@var{A}))} will likely be more
91 ## @seealso{svd, eigs}
93 function [u, s, v, flag] = svds (A, k, sigma, opts)
95 persistent root2 = sqrt (2);
97 if (nargin < 1 || nargin > 4)
102 error ("svds: A must be a 2D matrix");
106 opts.tol = 1e-10 / root2;
110 if (!isstruct (opts))
111 error ("svds: OPTS must be a structure");
113 if (!isfield (opts, "tol"))
114 opts.tol = 1e-10 / root2;
116 opts.tol = opts.tol / root2;
118 if (isfield (opts, "v0"))
119 if (!isvector (opts.v0) || (length (opts.v0) != sum (size (A))))
120 error ("svds: OPTS.v0 must be a vector with rows(A)+columns(A) entries");
125 if (nargin < 3 || strcmp (sigma, "L"))
131 elseif (isscalar (sigma) && isnumeric (sigma) && isreal (sigma))
133 error ("svds: SIGMA must be a positive real value");
136 error ("svds: SIGMA must be a positive real value or the string 'L'");
140 max_a = max (abs (A(:)));
142 s = zeros (k, 1); # special case of zero matrix
150 ## Scale everything by the 1-norm to make things more stable.
153 ## Call to eigs is always a symmetric matrix by construction
155 b_opts.tol = opts.tol / max_a;
157 if (!ischar (b_sigma))
158 b_sigma = b_sigma / max_a;
162 ## Find the smallest eigenvalues
163 ## The eigenvalues returns by eigs for sigma=0 are symmetric about 0.
164 ## As we are only interested in the positive eigenvalues, we have to
165 ## double k and then throw out the k negative eigenvalues.
166 ## Separately, if sigma is non-zero, but smaller than the smallest
167 ## singular value, ARPACK may not return k eigenvalues. However, as
168 ## computation scales with k we'd like to avoid doubling k for all
169 ## scalar values of sigma.
172 b_k = k; # Normal case, find just the k largest eigenvalues
176 [V, s, flag] = eigs ([sparse(m,m), b; b', sparse(n,n)],
177 b_k, b_sigma, b_opts);
180 s = eigs ([sparse(m,m), b; b', sparse(n,n)], b_k, b_sigma, b_opts);
188 ## We wish to exclude all eigenvalues that are less than zero as these
189 ## are artifacts of the way the matrix passed to eigs is formed. There
190 ## is also the possibility that the value of sigma chosen is exactly
191 ## a singular value, and in that case we're dead!! So have to rely on
192 ## the warning from eigs. We exclude the singular values which are
193 ## less than or equal to zero to within some tolerance scaled by the
194 ## norm since if we don't we might end up with too many singular
196 tol = norma * opts.tol;
198 if (length (ind) < k)
199 ## Too few eigenvalues returned. Add in any zero eigenvalues of B,
200 ## including the nominally negative ones.
201 zind = find (abs (s) <= tol);
202 p = min (length (zind), k - length (ind));
203 ind = [ind; zind(1:p)];
204 elseif (length (ind) > k)
205 ## Too many eigenvalues returned. Select according to criterium.
207 ind = ind(end+1-k:end); # smallest eigenvalues
209 ind = ind(1:k); # largest eigenvalues
215 warning ("returning fewer singular values than requested");
217 warning ("try increasing the value of sigma");
232 u = root2 * V(1:m,ind);
234 v = root2 * V(m+1:end,ind);
238 flag = norm (A*v - u*s, 1) > root2 * opts.tol * norm (A, 1);
244 %!shared n, k, A, u, s, v, opts, rand_state, randn_state
247 %! A = sparse ([3:n,1:n,1:(n-2)],[1:(n-2),1:n,3:n],[ones(1,n-2),0.4*n*ones(1,n),ones(1,n-2)]);
248 %! [u,s,v] = svd (full (A));
250 %! [~, idx] = sort (abs(s));
254 %! randn_state = randn ("state");
255 %! rand_state = rand ("state");
256 %! randn ("state", 42); % Initialize to make normest function reproducible
257 %! rand ("state", 42);
258 %! opts.v0 = rand (2*n,1); % Initialize eigs ARPACK starting vector
259 %! % to guarantee reproducible results
262 %! [u2,s2,v2,flag] = svds (A,k);
264 %! assert (flag, !1);
265 %! assert (s2, s(end:-1:end-k+1), 1e-10);
267 %!testif HAVE_ARPACK, HAVE_UMFPACK
268 %! [u2,s2,v2,flag] = svds (A,k,0,opts);
270 %! assert (flag, !1);
271 %! assert (s2, s(k:-1:1), 1e-10);
273 %!testif HAVE_ARPACK, HAVE_UMFPACK
275 %! % Don't put sigma right on a singular value or there are convergence issues
276 %! sigma = 0.99*s(idx) + 0.01*s(idx+1);
277 %! [u2,s2,v2,flag] = svds (A,k,sigma,opts);
279 %! assert (flag, !1);
280 %! assert (s2, s((idx+floor(k/2)):-1:(idx-floor(k/2))), 1e-10);
283 %! [u2,s2,v2,flag] = svds (zeros (10), k);
284 %! assert (u2, eye (10, k));
285 %! assert (s2, zeros (k));
286 %! assert (v2, eye (10, 7));
289 %! s = svds (speye (10));
290 %! assert (s, ones (6, 1), 2*eps);
293 %! ## Restore random number generator seeds at end of tests
294 %! rand ("state", rand_state);
295 %! randn ("state", randn_state);