Golub-Kahan-Lanczos (SVDL)
Usage
#
IterativeSolvers.svdl — Function
svdl(A) -> Σ, L, [history]
Compute some singular values (and optionally vectors) using Golub-Kahan-Lanczos bidiagonalization [1] with thick restarting [2].
If log is set to true is given, method will output a tuple X, L, ch. Where ch is a ConvergenceHistory object. Otherwise it will only return X, L.
Arguments
-
A: The matrix or matrix-like object whose singular values are desired.
Keywords
-
nsv::Int = 6: number of singular values requested; -
v0 = random unit vector: starting guess vector in the domain ofA. The length ofqshould be the number of columns inA; -
k::Int = 2nsv: maximum number of Lanczos vectors to compute before restarting; -
j::Int = nsv: number of vectors to keep at the end of the restart. We don’t recommend j < nsv; -
maxiter::Int = minimum(size(A)): maximum number of iterations to run; -
verbose::Bool = false: print information at each iteration; -
tol::Real = √eps(): maximum absolute error in each desired singular value; -
reltol::Real=√eps(): maximum error in each desired singular value relative to the estimated norm of the input matrix; -
method::Symbol=:ritz: restarting algorithm to use. Valid choices are:-
:ritz: Thick restart with Ritz values [Wu2000]. -
:harmonic: Restart with harmonic Ritz values [Baglama2005].
-
-
vecs::Symbol = :none: singular vectors to return.-
:both: Both left and right singular vectors are returned. -
:left: Only the left singular vectors are returned. -
:right: Only the right singular vectors are returned. -
:none: No singular vectors are returned.
-
-
dolock::Bool=false: Iftrue, locks converged Ritz values, removing them from the Krylov subspace being searched in the next macroiteration; -
log::Bool = false: output an extra element of typeConvergenceHistorycontaining extra information of the method execution.
Return values
if log is false
-
Σ: list of the desired singular values ifvecs == :none(the default), otherwise returns anSVDobject with the desired singular vectors filled in; -
L: computed partial factorizations of A.
if log is true
-
Σ: list of the desired singular values ifvecs == :none(the default),
otherwise returns an SVD object with the desired singular vectors filled in;
-
L: computed partial factorizations of A; -
history: convergence history.
ConvergenceHistory keys
-
:betas=>betas: The history of the computed betas. -
:Bs=>Bs: The history of the computed projected matrices. -
:ritz=>ritzvalhist: Ritz values computed at each iteration. -
:conv=>convhist: Convergence data.
Implementation details
The implementation of thick restarting follows closely that of SLEPc as described in [3]. Thick restarting can be turned off by setting k = maxiter, but most of the time this is not desirable.
The singular vectors are computed directly by forming the Ritz vectors from the product of the Lanczos vectors L.P/L.Q and the singular vectors of L.B. Additional accuracy in the singular triples can be obtained using inverse iteration.