BiLQ
Krylov.bilq
— Function(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 the BiLQ method.
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 the method
(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.
Reference
- A. Montoison and D. Orban, BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property, SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145–1166, 2020.
Krylov.bilq!
— Functionsolver = 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
.
QMR
Krylov.qmr
— Function(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 the QMR method.
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 the method
(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
- R. W. Freund and N. M. Nachtigal, QMR : a quasi-minimal residual method for non-Hermitian linear systems, Numerische mathematik, Vol. 60(1), pp. 315–339, 1991.
- R. W. Freund and N. M. Nachtigal, An implementation of the QMR method based on coupled two-term recurrences, SIAM Journal on Scientific Computing, Vol. 15(2), pp. 313–337, 1994.
- A. Montoison and D. Orban, BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property, SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145–1166, 2020.
Krylov.qmr!
— Functionsolver = 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
.
USYMLQ
Krylov.usymlq
— Function(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 the method
(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
- M. A. Saunders, H. D. Simon, and E. L. Yip, Two Conjugate-Gradient-Type Methods for Unsymmetric Linear Equations, SIAM Journal on Numerical Analysis, 25(4), pp. 927–940, 1988.
- A. Buttari, D. Orban, D. Ruiz and D. Titley-Peloquin, A tridiagonalization method for symmetric saddle-point and quasi-definite systems, SIAM Journal on Scientific Computing, 41(5), pp. 409–432, 2019.
- A. Montoison and D. Orban, BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property, SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145–1166, 2020.
Krylov.usymlq!
— Functionsolver = 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
.
USYMQR
Krylov.usymqr
— Function(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 the USYMQR method.
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 the method
(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
- M. A. Saunders, H. D. Simon, and E. L. Yip, Two Conjugate-Gradient-Type Methods for Unsymmetric Linear Equations, SIAM Journal on Numerical Analysis, 25(4), pp. 927–940, 1988.
- A. Buttari, D. Orban, D. Ruiz and D. Titley-Peloquin, A tridiagonalization method for symmetric saddle-point and quasi-definite systems, SIAM Journal on Scientific Computing, 41(5), pp. 409–432, 2019.
- A. Montoison and D. Orban, BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property, SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145–1166, 2020.
Krylov.usymqr!
— Functionsolver = 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
.
CGS
Krylov.cgs
— Function(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 the method
(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
- P. Sonneveld, CGS, A Fast Lanczos-Type Solver for Nonsymmetric Linear systems, SIAM Journal on Scientific and Statistical Computing, 10(1), pp. 36–52, 1989.
Krylov.cgs!
— Functionsolver = 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
.
BiCGSTAB
Krylov.bicgstab
— Function(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 the BICGSTAB method. 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 the method
(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
- H. A. van der Vorst, Bi-CGSTAB: A fast and smoothly converging variant of Bi-CG for the solution of nonsymmetric linear systems, SIAM Journal on Scientific and Statistical Computing, 13(2), pp. 631–644, 1992.
- G. L.G. Sleijpen and D. R. Fokkema, BiCGstab(ℓ) for linear equations involving unsymmetric matrices with complex spectrum, Electronic Transactions on Numerical Analysis, 1, pp. 11–32, 1993.
Krylov.bicgstab!
— Functionsolver = 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
.
DIOM
Krylov.diom
— Function(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 direct incomplete orthogonalization method.
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.
- Left preconditioning : M⁻¹Ax = M⁻¹b
- Right preconditioning : AN⁻¹u = b with x = N⁻¹u
- Split preconditioning : M⁻¹AN⁻¹u = M⁻¹b with x = N⁻¹u
DIOM can be warm-started from an initial guess x0
with the method
(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
- Y. Saad, Practical use of some krylov subspace methods for solving indefinite and nonsymmetric linear systems, SIAM journal on scientific and statistical computing, 5(1), pp. 203–228, 1984.
Krylov.diom!
— Functionsolver = 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
.
FOM
Krylov.fom
— Function(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 method.
FOM algorithm is based on the Arnoldi process and a Galerkin condition.
This implementation allows a left preconditioner M and a right preconditioner N.
- Left preconditioning : M⁻¹Ax = M⁻¹b
- Right preconditioning : AN⁻¹u = b with x = N⁻¹u
- Split preconditioning : M⁻¹AN⁻¹u = M⁻¹b with x = N⁻¹u
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 exceed memory
.
FOM can be warm-started from an initial guess x0
with the method
(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
- Y. Saad, Krylov subspace methods for solving unsymmetric linear systems, Mathematics of computation, Vol. 37(155), pp. 105–126, 1981.
Krylov.fom!
— Functionsolver = 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
.
DQGMRES
Krylov.dqgmres
— Function(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 method.
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.
- Left preconditioning : M⁻¹Ax = M⁻¹b
- Right preconditioning : AN⁻¹u = b with x = N⁻¹u
- Split preconditioning : M⁻¹AN⁻¹u = M⁻¹b with x = N⁻¹u
DQGMRES can be warm-started from an initial guess x0
with the method
(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
- Y. Saad and K. Wu, DQGMRES: a quasi minimal residual algorithm based on incomplete orthogonalization, Numerical Linear Algebra with Applications, Vol. 3(4), pp. 329–343, 1996.
Krylov.dqgmres!
— Functionsolver = 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
.
GMRES
Krylov.gmres
— Function(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 method.
GMRES algorithm is based on the Arnoldi process and computes a sequence of approximate solutions with the minimal residual property.
This implementation allows a left preconditioner M and a right preconditioner N.
- Left preconditioning : M⁻¹Ax = M⁻¹b
- Right preconditioning : AN⁻¹u = b with x = N⁻¹u
- Split preconditioning : M⁻¹AN⁻¹u = M⁻¹b with x = N⁻¹u
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 exceed memory
.
GMRES can be warm-started from an initial guess x0
with the method
(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
- Y. Saad and M. H. Schultz, GMRES: A Generalized Minimal Residual Algorithm for Solving Nonsymmetric Linear Systems, SIAM Journal on Scientific and Statistical Computing, Vol. 7(3), pp. 856–869, 1986.
Krylov.gmres!
— Functionsolver = 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
.