CG

Krylov.cgFunction
(x, stats) = cg(A, b::AbstractVector{FC};
                M=I, atol::T=√eps(T), rtol::T=√eps(T),
                itmax::Int=0, radius::T=zero(T), linesearch::Bool=false,
                verbose::Int=0, history::Bool=false,
                ldiv::Bool=false, callback=solver->false)

T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.

The conjugate gradient method to solve the symmetric linear system Ax = b.

The method does not abort if A is not definite.

A preconditioner M may be provided in the form of a linear operator and is assumed to be symmetric and positive definite. M also indicates the weighted norm in which residuals are measured.

If itmax=0, the default number of iterations is set to 2 * n, with n = length(b).

CG can be warm-started from an initial guess x0 with

(x, stats) = cg(A, b, x0; kwargs...)

where kwargs are the same keyword arguments as above.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

Reference

source
Krylov.cg!Function
solver = cg!(solver::CgSolver, A, b; kwargs...)
solver = cg!(solver::CgSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of cg.

See CgSolver for more details about the solver.

source

CR

Krylov.crFunction
(x, stats) = cr(A, b::AbstractVector{FC};
                M=I, atol::T=√eps(T), rtol::T=√eps(T), γ::T=√eps(T), itmax::Int=0,
                radius::T=zero(T), verbose::Int=0, linesearch::Bool=false, history::Bool=false,
                ldiv::Bool=false, callback=solver->false)

T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.

A truncated version of Stiefel’s Conjugate Residual method to solve the symmetric linear system Ax = b or the least-squares problem min ‖b - Ax‖. The matrix A must be positive semi-definite.

A preconditioner M may be provided in the form of a linear operator and is assumed to be symmetric and positive definite. M also indicates the weighted norm in which residuals are measured.

In a linesearch context, 'linesearch' must be set to 'true'.

If itmax=0, the default number of iterations is set to 2 * n, with n = length(b).

CR can be warm-started from an initial guess x0 with

(x, stats) = cr(A, b, x0; kwargs...)

where kwargs are the same keyword arguments as above.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

References

source
Krylov.cr!Function
solver = cr!(solver::CrSolver, A, b; kwargs...)
solver = cr!(solver::CrSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of cr.

See CrSolver for more details about the solver.

source

CG-LANCZOS

Krylov.cg_lanczosFunction
(x, stats) = cg_lanczos(A, b::AbstractVector{FC};
                        M=I, atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
                        check_curvature::Bool=false, verbose::Int=0, history::Bool=false,
                        ldiv::Bool=false, callback=solver->false)

T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.

The Lanczos version of the conjugate gradient method to solve the symmetric linear system

Ax = b

The method does not abort if A is not definite.

A preconditioner M may be provided in the form of a linear operator and is assumed to be hermitian and positive definite.

CG-LANCZOS can be warm-started from an initial guess x0 with

(x, stats) = cg_lanczos(A, b, x0; kwargs...)

where kwargs are the same keyword arguments as above.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

References

source
Krylov.cg_lanczos!Function
solver = cg_lanczos!(solver::CgLanczosSolver, A, b; kwargs...)
solver = cg_lanczos!(solver::CgLanczosSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of cg_lanczos.

See CgLanczosSolver for more details about the solver.

source

CG-LANCZOS-SHIFT

Krylov.cg_lanczos_shiftFunction
(x, stats) = cg_lanczos_shift(A, b::AbstractVector{FC}, shifts::AbstractVector{T};
                              M=I, atol::T=√eps(T), rtol::T=√eps(T),
                              itmax::Int=0, check_curvature::Bool=false,
                              verbose::Int=0, history::Bool=false,
                              ldiv::Bool=false, callback=solver->false)

T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.

The Lanczos version of the conjugate gradient method to solve a family of shifted systems

(A + αI) x = b  (α = α₁, ..., αₙ)

The method does not abort if A + αI is not definite.

A preconditioner M may be provided in the form of a linear operator and is assumed to be hermitian and positive definite.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

source