BiLQ

Krylov.bilqFunction
(x, stats) = bilq(A, b::AbstractVector{FC}; c::AbstractVector{FC}=b,
                  atol::T=√eps(T), rtol::T=√eps(T), transfer_to_bicg::Bool=true,
                  itmax::Int=0, verbose::Int=0, history::Bool=false,
                  callback=solver->false)

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

Solve the square linear system Ax = b using BiLQ.

BiLQ is based on the Lanczos biorthogonalization process and requires two initial vectors b and c. The relation bᴴc ≠ 0 must be satisfied and by default c = b. When A is symmetric and b = c, BiLQ is equivalent to SYMMLQ.

An option gives the possibility of transferring to the BiCG point, when it exists. The transfer is based on the residual norm.

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

(x, stats) = bilq(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.bilq!Function
solver = bilq!(solver::BilqSolver, A, b; kwargs...)
solver = bilq!(solver::BilqSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of bilq.

See BilqSolver for more details about the solver.

source

QMR

Krylov.qmrFunction
(x, stats) = qmr(A, b::AbstractVector{FC}; c::AbstractVector{FC}=b,
                 atol::T=√eps(T), rtol::T=√eps(T),
                 itmax::Int=0, verbose::Int=0, history::Bool=false,
                 callback=solver->false)

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

Solve the square linear system Ax = b using QMR.

QMR is based on the Lanczos biorthogonalization process and requires two initial vectors b and c. The relation bᴴc ≠ 0 must be satisfied and by default c = b. When A is symmetric and b = c, QMR is equivalent to MINRES.

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

(x, stats) = qmr(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.qmr!Function
solver = qmr!(solver::QmrSolver, A, b; kwargs...)
solver = qmr!(solver::QmrSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of qmr.

See QmrSolver for more details about the solver.

source

USYMLQ

Krylov.usymlqFunction
(x, stats) = usymlq(A, b::AbstractVector{FC}, c::AbstractVector{FC};
                    atol::T=√eps(T), rtol::T=√eps(T), transfer_to_usymcg::Bool=true,
                    itmax::Int=0, verbose::Int=0, history::Bool=false,
                    callback=solver->false)

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

Solve the linear system Ax = b using the USYMLQ method.

USYMLQ is based on the orthogonal tridiagonalization process and requires two initial nonzero vectors b and c. The vector c is only used to initialize the process and a default value can be b or Aᴴb depending on the shape of A. The error norm ‖x - x*‖ monotonously decreases in USYMLQ. It's considered as a generalization of SYMMLQ.

It can also be applied to under-determined and over-determined problems. In all cases, problems must be consistent.

An option gives the possibility of transferring to the USYMCG point, when it exists. The transfer is based on the residual norm.

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

(x, stats) = usymlq(A, b, c, 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.usymlq!Function
solver = usymlq!(solver::UsymlqSolver, A, b, c; kwargs...)
solver = usymlq!(solver::UsymlqSolver, A, b, c, x0; kwargs...)

where kwargs are keyword arguments of usymlq.

See UsymlqSolver for more details about the solver.

source

USYMQR

Krylov.usymqrFunction
(x, stats) = usymqr(A, b::AbstractVector{FC}, c::AbstractVector{FC};
                    atol::T=√eps(T), rtol::T=√eps(T),
                    itmax::Int=0, verbose::Int=0, history::Bool=false,
                    callback=solver->false)

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

Solve the linear system Ax = b using USYMQR.

USYMQR is based on the orthogonal tridiagonalization process and requires two initial nonzero vectors b and c. The vector c is only used to initialize the process and a default value can be b or Aᴴb depending on the shape of A. The residual norm ‖b - Ax‖ monotonously decreases in USYMQR. It's considered as a generalization of MINRES.

It can also be applied to under-determined and over-determined problems. USYMQR finds the minimum-norm solution if problems are inconsistent.

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

(x, stats) = usymqr(A, b, c, 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.usymqr!Function
solver = usymqr!(solver::UsymqrSolver, A, b, c; kwargs...)
solver = usymqr!(solver::UsymqrSolver, A, b, c, x0; kwargs...)

where kwargs are keyword arguments of usymqr.

See UsymqrSolver for more details about the solver.

source

CGS

Krylov.cgsFunction
(x, stats) = cgs(A, b::AbstractVector{FC}; c::AbstractVector{FC}=b,
                 M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                 itmax::Int=0, 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}.

Solve the consistent linear system Ax = b using conjugate gradient squared algorithm. CGS requires two initial vectors b and c. The relation bᴴc ≠ 0 must be satisfied and by default c = b.

From "Iterative Methods for Sparse Linear Systems (Y. Saad)" :

«The method is based on a polynomial variant of the conjugate gradients algorithm. Although related to the so-called bi-conjugate gradients (BCG) algorithm, it does not involve adjoint matrix-vector multiplications, and the expected convergence rate is about twice that of the BCG algorithm.

The Conjugate Gradient Squared algorithm works quite well in many cases. However, one difficulty is that, since the polynomials are squared, rounding errors tend to be more damaging than in the standard BCG algorithm. In particular, very high variations of the residual vectors often cause the residual norms computed to become inaccurate.

TFQMR and BICGSTAB were developed to remedy this difficulty.»

This implementation allows a left preconditioner M and a right preconditioner N.

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

(x, stats) = cgs(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.cgs!Function
solver = cgs!(solver::CgsSolver, A, b; kwargs...)
solver = cgs!(solver::CgsSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of cgs.

See CgsSolver for more details about the solver.

source

BiCGSTAB

Krylov.bicgstabFunction
(x, stats) = bicgstab(A, b::AbstractVector{FC}; c::AbstractVector{FC}=b,
                      M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                      itmax::Int=0, 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}.

Solve the square linear system Ax = b using BICGSTAB. BICGSTAB requires two initial vectors b and c. The relation bᴴc ≠ 0 must be satisfied and by default c = b.

The Biconjugate Gradient Stabilized method is a variant of BiCG, like CGS, but using different updates for the Aᴴ-sequence in order to obtain smoother convergence than CGS.

If BICGSTAB stagnates, we recommend DQGMRES and BiLQ as alternative methods for unsymmetric square systems.

BICGSTAB stops when itmax iterations are reached or when ‖rₖ‖ ≤ atol + ‖b‖ * rtol. atol is an absolute tolerance and rtol is a relative tolerance.

Additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every verbose iterations.

This implementation allows a left preconditioner M and a right preconditioner N.

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

(x, stats) = bicgstab(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.bicgstab!Function
solver = bicgstab!(solver::BicgstabSolver, A, b; kwargs...)
solver = bicgstab!(solver::BicgstabSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of bicgstab.

See BicgstabSolver for more details about the solver.

source

DIOM

Krylov.diomFunction
(x, stats) = diom(A, b::AbstractVector{FC}; memory::Int=20,
                  M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                  reorthogonalization::Bool=false, itmax::Int=0,
                  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}.

Solve the consistent linear system Ax = b using DIOM.

DIOM only orthogonalizes the new vectors of the Krylov basis against the memory most recent vectors. If CG is well defined on Ax = b and memory = 2, DIOM is theoretically equivalent to CG. If k ≤ memory where k is the number of iterations, DIOM is theoretically equivalent to FOM. Otherwise, DIOM interpolates between CG and FOM and is similar to CG with partial reorthogonalization.

Partial reorthogonalization is available with the reorthogonalization option.

An advantage of DIOM is that nonsymmetric or symmetric indefinite or both nonsymmetric and indefinite systems of linear equations can be handled by this single algorithm.

This implementation allows a left preconditioner M and a right preconditioner N.

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

(x, stats) = diom(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.diom!Function
solver = diom!(solver::DiomSolver, A, b; kwargs...)
solver = diom!(solver::DiomSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of diom.

Note that the memory keyword argument is the only exception. It's required to create a DiomSolver and can't be changed later.

See DiomSolver for more details about the solver.

source

FOM

Krylov.fomFunction
(x, stats) = fom(A, b::AbstractVector{FC}; memory::Int=20,
                 M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                 reorthogonalization::Bool=false, itmax::Int=0,
                 restart::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}.

Solve the linear system Ax = b using FOM.

FOM algorithm is based on the Arnoldi process and a Galerkin condition.

This implementation allows a left preconditioner M and a right preconditioner N. Full reorthogonalization is available with the reorthogonalization option.

If restart = true, the restarted version FOM(k) is used with k = memory. If restart = false, the parameter memory should be used as a hint of the number of iterations to limit dynamic memory allocations. More storage will be allocated only if the number of iterations exceeds memory.

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

(x, stats) = fom(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.fom!Function
solver = fom!(solver::FomSolver, A, b; kwargs...)
solver = fom!(solver::FomSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of fom.

Note that the memory keyword argument is the only exception. It's required to create a FomSolver and can't be changed later.

See FomSolver for more details about the solver.

source

DQGMRES

Krylov.dqgmresFunction
(x, stats) = dqgmres(A, b::AbstractVector{FC}; memory::Int=20,
                     M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                     reorthogonalization::Bool=false, itmax::Int=0,
                     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}.

Solve the consistent linear system Ax = b using DQGMRES.

DQGMRES algorithm is based on the incomplete Arnoldi orthogonalization process and computes a sequence of approximate solutions with the quasi-minimal residual property.

DQGMRES only orthogonalizes the new vectors of the Krylov basis against the memory most recent vectors. If MINRES is well defined on Ax = b and memory = 2, DQGMRES is theoretically equivalent to MINRES. If k ≤ memory where k is the number of iterations, DQGMRES is theoretically equivalent to GMRES. Otherwise, DQGMRES interpolates between MINRES and GMRES and is similar to MINRES with partial reorthogonalization.

Partial reorthogonalization is available with the reorthogonalization option.

This implementation allows a left preconditioner M and a right preconditioner N.

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

(x, stats) = dqgmres(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.dqgmres!Function
solver = dqgmres!(solver::DqgmresSolver, A, b; kwargs...)
solver = dqgmres!(solver::DqgmresSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of dqgmres.

Note that the memory keyword argument is the only exception. It's required to create a DqgmresSolver and can't be changed later.

See DqgmresSolver for more details about the solver.

source

GMRES

Krylov.gmresFunction
(x, stats) = gmres(A, b::AbstractVector{FC}; memory::Int=20,
                   M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                   reorthogonalization::Bool=false, itmax::Int=0,
                   restart::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}.

Solve the linear system Ax = b using GMRES.

GMRES algorithm is based on the Arnoldi process and computes a sequence of approximate solutions with the minimum residual.

This implementation allows a left preconditioner M and a right preconditioner N. Full reorthogonalization is available with the reorthogonalization option.

If restart = true, the restarted version GMRES(k) is used with k = memory. If restart = false, the parameter memory should be used as a hint of the number of iterations to limit dynamic memory allocations. More storage will be allocated only if the number of iterations exceeds memory.

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

(x, stats) = gmres(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.gmres!Function
solver = gmres!(solver::GmresSolver, A, b; kwargs...)
solver = gmres!(solver::GmresSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of gmres.

Note that the memory keyword argument is the only exception. It's required to create a GmresSolver and can't be changed later.

See GmresSolver for more details about the solver.

source

FGMRES

Krylov.fgmresFunction
(x, stats) = fgmres(A, b::AbstractVector{FC}; memory::Int=20,
                    M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
                    reorthogonalization::Bool=false, itmax::Int=0,
                    restart::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}.

Solve the linear system Ax = b using FGMRES.

FGMRES computes a sequence of approximate solutions with minimum residual. FGMRES is a variant of GMRES that allows changes in the right preconditioner at each iteration.

This implementation allows a left preconditioner M and a flexible right preconditioner N. A situation in which the preconditioner is "not constant" is when a relaxation-type method, a Chebyshev iteration or another Krylov subspace method is used as a preconditioner. Compared to GMRES, there is no additional cost incurred in the arithmetic but the memory requirement almost doubles. Thus, GMRES is recommended if the right preconditioner N is constant.

Full reorthogonalization is available with the reorthogonalization option.

If restart = true, the restarted version FGMRES(k) is used with k = memory. If restart = false, the parameter memory should be used as a hint of the number of iterations to limit dynamic memory allocations. More storage will be allocated only if the number of iterations exceeds memory.

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

(x, stats) = fgmres(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.fgmres!Function
solver = fgmres!(solver::FgmresSolver, A, b; kwargs...)
solver = fgmres!(solver::FgmresSolver, A, b, x0; kwargs...)

where kwargs are keyword arguments of fgmres.

Note that the memory keyword argument is the only exception. It's required to create a FgmresSolver and can't be changed later.

See FgmresSolver for more details about the solver.

source