Available algorithms
Orthogonalization algorithms
KrylovKit.ClassicalGramSchmidt
— Type.ClassicalGramSchmidt()
Represents the classical Gram Schmidt algorithm for orthogonalizing different vectors, typically not an optimal choice.
KrylovKit.ModifiedGramSchmidt
— Type.ModifiedGramSchmidt()
Represents the modified Gram Schmidt algorithm for orthogonalizing different vectors, typically a reasonable choice for linear systems but not for eigenvalue solvers with a large Krylov dimension.
KrylovKit.ClassicalGramSchmidt2
— Type.ClassicalGramSchmidt2()
Represents the classical Gram Schmidt algorithm with a second reorthogonalization step always taking place.
KrylovKit.ModifiedGramSchmidt2
— Type.ModifiedGramSchmidt2()
Represents the modified Gram Schmidt algorithm with a second reorthogonalization step always taking place.
KrylovKit.ClassicalGramSchmidtIR
— Type.ClassicalGramSchmidtIR(η::Real)
Represents the classical Gram Schmidt algorithm with zero or more reorthogonalization steps being applied untill the norm of the vector after an orthogonalization step has not decreased by a factor smaller than η
with respect to the norm before the step.
KrylovKit.ModifiedGramSchmidtIR
— Type.ModifiedGramSchmidtIR(η::Real)
Represents the modified Gram Schmidt algorithm with zero or more reorthogonalization steps being applied untill the norm of the vector after an orthogonalization step has not decreased by a factor smaller than η
with respect to the norm before the step.
General Krylov algorithms
KrylovKit.Lanczos
— Type.Lanczos(orth::Orthogonalizer = KrylovDefaults.orth; krylovdim = KrylovDefaults.krylovdim,
maxiter::Int = KrylovDefaults.maxiter, tol = KrylovDefaults.tol)
Represents the Lanczos algorithm for building the Krylov subspace; assumes the linear operator is real symmetric or complex Hermitian. Can be used in eigsolve
and exponentiate
. The corresponding algorithms will build a Krylov subspace of size at most krylovdim
, which will be repeated at most maxiter
times and will stop when the norm of the residual of the Lanczos factorization is smaller than tol
. The orthogonalizer orth
will be used to orthogonalize the different Krylov vectors.
Use Arnoldi
for non-symmetric or non-Hermitian linear operators.
See also: factorize
, eigsolve
, exponentiate
, Arnoldi
, Orthogonalizer
KrylovKit.Arnoldi
— Type.Arnoldi(orth::Orthogonalizer = KrylovDefaults.orth; krylovdim = KrylovDefaults.krylovdim,
maxiter::Int = KrylovDefaults.maxiter, tol = KrylovDefaults.tol)
Represents the Arnoldi algorithm for building the Krylov subspace for a general matrix or linear operator. Can be used in eigsolve
and exponentiate
. The corresponding algorithms will build a Krylov subspace of size at most krylovdim
, which will be repeated at most maxiter
times and will stop when the norm of the residual of the Arnoldi factorization is smaller than tol
. The orthogonalizer orth
will be used to orthogonalize the different Krylov vectors.
Use Lanczos
for real symmetric or complex Hermitian linear operators.
See also: eigsolve
, exponentiate
, Lanczos
, Orthogonalizer
Specific algorithms for linear systems
KrylovKit.CG
— Type.CG(; maxiter = KrylovDefaults.maxiter, atol = 0, rtol = KrylovDefaults.tol)
Construct an instance of the conjugate gradient algorithm with specified parameters, which can be passed to linsolve
in order to iteratively solve a linear system with a positive definite (and thus symmetric or hermitian) coefficent matrix or operator. The CG
method will search for the optimal x
in a Krylov subspace of maximal size maxiter
, or stop when norm(A*x - b) < max(atol, rtol*norm(b))
.
KrylovKit.MINRES
— Type.MINRES(; maxiter = KrylovDefaults.maxiter, atol = 0, rtol = KrylovDefaults.tol)
Construct an instance of the conjugate gradient algorithm with specified parameters, which can be passed to linsolve
in order to iteratively solve a linear system with a real symmetric or complex hermitian coefficent matrix or operator. The MINRES
method will search for the optimal x
in a Krylov subspace of maximal size maxiter
, or stop when norm(A*x - b) < max(atol, rtol*norm(b))
.
KrylovKit.GMRES
— Type.GMRES(orth::Orthogonalizer = KrylovDefaults.orth; maxiter = KrylovDefaults.maxiter,
krylovdim = KrylovDefaults.krylovdim, atol = 0, rtol = KrylovDefaults.tol)
Construct an instance of the GMRES algorithm with specified parameters, which can be passed to linsolve
in order to iteratively solve a linear system. The GMRES
method will search for the optimal x
in a Krylov subspace of maximal size krylovdim
, and repeat this process for at most maxiter
times, or stop when norm(A*x - b) < max(atol, rtol*norm(b))
.
In building the Krylov subspace, GMRES
will use the orthogonalizer orth
.
Note that we do not follow the nomenclature in the traditional literature on GMRES
, where krylovdim
is referred to as the restart parameter, and every new Krylov vector counts as an iteration. I.e. our iteration count should rougly be multiplied by krylovdim
to obtain the conventional iteration count.
KrylovKit.BiCG
— Type.BiCG(; maxiter = KrylovDefaults.maxiter, atol = 0, rtol = KrylovDefaults.tol)
Construct an instance of the Biconjugate gradient algorithm with specified parameters, which can be passed to linsolve
in order to iteratively solve a linear system general linear map, of which the adjoint can also be applied. The BiCG
method will search for the optimal x
in a Krylov subspace of maximal size maxiter
, or stop when norm(A*x - b) < max(atol, rtol*norm(b))
.
KrylovKit.BiCGStab
— Type.BiCGStab(; maxiter = KrylovDefaults.maxiter, atol = 0, rtol = KrylovDefaults.tol)
Construct an instance of the Biconjugate gradient algorithm with specified parameters, which can be passed to linsolve
in order to iteratively solve a linear system general linear map. The BiCGStab
method will search for the optimal x
in a Krylov subspace of maximal size maxiter
, or stop when norm(A*x - b) < max(atol, rtol*norm(b))
.
Specific algorithms for sigular values
KrylovKit.GKL
— Type.GKL(orth::Orthogonalizer = KrylovDefaults.orth; krylovdim = KrylovDefaults.krylovdim,
maxiter::Int = KrylovDefaults.maxiter, tol = KrylovDefaults.tol)
Represents the Golub-Kahan-Lanczos bidiagonalization algorithm for sequentially building a Krylov-like factorization of a genereal matrix or linear operator with a bidiagonal reduced matrix. Can be used in svdsolve
. The corresponding algorithm builds a Krylov subspace of size at most krylovdim
, which will be repeated at most maxiter
times and will stop when the norm of the residual of the Arnoldi factorization is smaller than tol
. The orthogonalizer orth
will be used to orthogonalize the different Krylov vectors.
See also: svdsolve
, Orthogonalizer
Default values
KrylovKit.KrylovDefaults
— Module.module KrylovDefaults
const orth = KrylovKit.ModifiedGramSchmidtIR(1/sqrt(2))
const krylovdim = 30
const maxiter = 100
const tol = 1e-12
end
A module listing the default values for the typical parameters in Krylov based algorithms:
orth
: the orthogonalization routine used to orthogonalize the Krylov basis in theLanczos
orArnoldi
iterationkrylovdim
: the maximal dimension of the Krylov subspace that will be constructedmaxiter
: the maximal number of outer iterations, i.e. the maximum number of times the Krylov subspace may be rebuilttol
: the tolerance to which the problem must be solved, based on a suitable error measure, e.g. the norm of some residual. !!! warning The default value oftol
is aFloat64
value, if you solve problems inFloat32
orComplexF32
arithmetic, you should always specify a newtol
as the default value will not be attainable.