Available algorithms

Orthogonalization algorithms

KrylovKit.ClassicalGramSchmidtType
ClassicalGramSchmidt()

Represents the classical Gram Schmidt algorithm for orthogonalizing different vectors, typically not an optimal choice.

source
KrylovKit.ModifiedGramSchmidtType
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.

source
KrylovKit.ClassicalGramSchmidtIRType
ClassicalGramSchmidtIR(η::Real = 1/sqrt(2))

Represents the classical Gram Schmidt algorithm with iterative (i.e. zero or more) reorthogonalization until 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. The default value corresponds to the Daniel-Gragg-Kaufman-Stewart condition.

source
KrylovKit.ModifiedGramSchmidtIRType
ModifiedGramSchmidtIR(η::Real = 1/sqrt(2))

Represents the modified Gram Schmidt algorithm with iterative (i.e. zero or more) reorthogonalization until 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. The default value corresponds to the Daniel-Gragg-Kaufman-Stewart condition.

source

General Krylov algorithms

KrylovKit.LanczosType
Lanczos(; krylovdim = KrylovDefaults.krylovdim, maxiter = KrylovDefaults.maxiter,
    tol = KrylovDefaults.tol, orth = KrylovDefaults.orth, eager = false, verbosity = 0)

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. Eager mode, as selected by eager = true, means that the algorithm that uses this Lanczos process (e.g. eigsolve) can try to finish its computation before the total Krylov subspace of dimension krylovdim is constructed. Default verbosity level verbosity is zero, meaning that no output will be printed.

Use Arnoldi for non-symmetric or non-Hermitian linear operators.

See also: factorize, eigsolve, exponentiate, Arnoldi, Orthogonalizer

source
KrylovKit.ArnoldiType
Arnoldi(; krylovdim = KrylovDefaults.krylovdim, maxiter = KrylovDefaults.maxiter,
    tol = KrylovDefaults.tol, orth = KrylovDefaults.orth, eager = false, verbosity = 0)

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. Eager mode, as selected by eager = true, means that the algorithm that uses this Arnoldi process (e.g. eigsolve) can try to finish its computation before the total Krylov subspace of dimension krylovdim is constructed. Default verbosity level verbosity is zero, meaning that no output will be printed.

Use Lanczos for real symmetric or complex Hermitian linear operators.

See also: eigsolve, exponentiate, Lanczos, Orthogonalizer

source

Specific algorithms for linear systems

KrylovKit.CGType
CG(; maxiter = KrylovDefaults.maxiter, tol = 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) coefficient 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) < tol. Default verbosity level verbosity is zero, meaning that no output will be printed.

See also: linsolve, MINRES, GMRES, BiCG, BiCGStab

source
KrylovKit.MINRESType
MINRES(; maxiter = KrylovDefaults.maxiter, tol = KrylovDefaults.tol)

!!! warning "Not implemented yet"

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 coefficient 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) < tol`. In building the Krylov subspace, `MINRES` will use the
orthogonalizer `orth`. Default verbosity level `verbosity` is zero, meaning that no
output will be printed.

See also: linsolve, CG, GMRES, BiCG, BiCGStab

source
KrylovKit.GMRESType
GMRES(; krylovdim = KrylovDefaults.krylovdim, maxiter = KrylovDefaults.maxiter,
    tol = KrylovDefaults.tol, orth::Orthogonalizer = KrylovDefaults.orth)

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) < tol. In building the Krylov subspace, GMRES will use the orthogonalizer orth. Default verbosity level verbosity is zero, meaning that no output will be printed.

Note that in the traditional nomenclature of GMRES, the parameter krylovdim is referred to as the restart parameter, and maxiter is the number of outer iterations, i.e. restart cycles. The total iteration count, i.e. the number of expansion steps, is roughly krylovdim times the number of iterations.

See also: linsolve, BiCG, BiCGStab, CG, MINRES

source
KrylovKit.BiCGType
BiCG(; maxiter = KrylovDefaults.maxiter, tol = KrylovDefaults.tol)

!!! warning "Not implemented yet"

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) < tol`. Default verbosity level `verbosity` is zero, meaning that no output will be
printed.

See also: linsolve, GMRES, CG, BiCGStab, MINRES

source
KrylovKit.BiCGStabType
BiCGStab(; maxiter = KrylovDefaults.maxiter, tol = 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) < tol`. Default verbosity level
`verbosity` is zero, meaning that no output will be printed.

See also: linsolve, GMRES, CG, BiCG, MINRES

source

Specific algorithms for generalized eigenvalue problems

KrylovKit.GolubYeType
GolubYe(; krylovdim = KrylovDefaults.krylovdim, maxiter = KrylovDefaults.maxiter,
    tol = KrylovDefaults.tol, orth = KrylovDefaults.orth, verbosity = 0)

Represents the Golub-Ye algorithm for solving hermitian (symmetric) generalized eigenvalue problems A x = λ B x with positive definite B, without the need for inverting B. Builds a Krylov subspace of size krylovdim starting from an estimate x by acting with (A - ρ(x) B), where ρ(x) = dot(x, A*x)/dot(x, B*x), and employing the Lanczos algorithm. This process is repeated at most maxiter times. In every iteration k>1, the subspace will also be expanded to size krylovdim+1 by adding $x_k - x_{k-1}$, which is known as the LOPCG correction and was suggested by Money and Ye. With krylovdim = 2, this algorithm becomes equivalent to LOPCG.

source

Specific algorithms for singular value problems

KrylovKit.GKLType
GKL(; krylovdim = KrylovDefaults.krylovdim, maxiter = KrylovDefaults.maxiter,
    tol = KrylovDefaults.tol, orth = KrylovDefaults.orth, verbosity = 0)

Represents the Golub-Kahan-Lanczos bidiagonalization algorithm for sequentially building a Krylov-like factorization of a general 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. Default verbosity level verbosity is zero, meaning that no output will be printed.

See also: svdsolve, Orthogonalizer

source

Default values

KrylovKit.KrylovDefaultsModule
module KrylovDefaults
    const orth = KrylovKit.ModifiedGramSchmidtIR()
    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 = ModifiedGramSchmidtIR(): the orthogonalization routine used to orthogonalize the Krylov basis in the Lanczos or Arnoldi iteration
  • krylovdim = 30: the maximal dimension of the Krylov subspace that will be constructed
  • maxiter = 100: the maximal number of outer iterations, i.e. the maximum number of times the Krylov subspace may be rebuilt
  • tol = 1e-12: 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 of tol is a Float64 value, if you solve problems in Float32 or ComplexF32 arithmetic, you should always specify a new tol as the default value will not be attainable.

source