BiLQ
Krylov.bilq — Function(x, stats) = bilq(A, b::AbstractVector{FC};
c::AbstractVector{FC}=b, transfer_to_bicg::Bool=true,
M=I, N=I, ldiv::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0, timemax::Float64=Inf,
verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = bilq(A, b, x0::AbstractVector; kwargs...)BiLQ can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the square linear system Ax = b of size n 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 Hermitian and b = c, BiLQ is equivalent to SYMMLQ. BiLQ requires support for adjoint(M) and adjoint(N) if preconditioners are provided.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :bilq.
For an in-place variant that reuses memory across solves, see bilq!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
c: the second initial vector of lengthnrequired by the Lanczos biorthogonalization process;transfer_to_bicg: transfer from the BiLQ point to the BiCG point, when it exists. The transfer is based on the residual norm;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
References
- 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.
- R. Fletcher, Conjugate gradient methods for indefinite systems, Numerical Analysis, Springer, pp. 73–89, 1976.
Krylov.bilq! — Functionworkspace = bilq!(workspace::BilqWorkspace, A, b; kwargs...)
workspace = bilq!(workspace::BilqWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of bilq.
See BilqWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :bilq to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
QMR
Krylov.qmr — Function(x, stats) = qmr(A, b::AbstractVector{FC};
c::AbstractVector{FC}=b, M=I, N=I, ldiv::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0, timemax::Float64=Inf, verbose::Int=0,
history::Bool=false, callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = qmr(A, b, x0::AbstractVector; kwargs...)QMR can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the square linear system Ax = b of size n 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 Hermitian and b = c, QMR is equivalent to MINRES. QMR requires support for adjoint(M) and adjoint(N) if preconditioners are provided.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :qmr.
For an in-place variant that reuses memory across solves, see qmr!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
c: the second initial vector of lengthnrequired by the Lanczos biorthogonalization process;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = qmr!(workspace::QmrWorkspace, A, b; kwargs...)
workspace = qmr!(workspace::QmrWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of qmr.
See QmrWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :qmr to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
CGS
Krylov.cgs — Function(x, stats) = cgs(A, b::AbstractVector{FC};
c::AbstractVector{FC}=b, M=I, N=I,
ldiv::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = cgs(A, b, x0::AbstractVector; kwargs...)CGS can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the consistent linear system Ax = b of size n using CGS. 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.»
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :cgs.
For an in-place variant that reuses memory across solves, see cgs!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
c: the second initial vector of lengthnrequired by the Lanczos biorthogonalization process;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = cgs!(workspace::CgsWorkspace, A, b; kwargs...)
workspace = cgs!(workspace::CgsWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of cgs.
See CgsWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :cgs to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
BiCGSTAB
Krylov.bicgstab — Function(x, stats) = bicgstab(A, b::AbstractVector{FC};
c::AbstractVector{FC}=b, M=I, N=I,
ldiv::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = bicgstab(A, b, x0::AbstractVector; kwargs...)BICGSTAB can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the square linear system Ax = b of size n 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :bicgstab.
For an in-place variant that reuses memory across solves, see bicgstab!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
c: the second initial vector of lengthnrequired by the Lanczos biorthogonalization process;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = bicgstab!(workspace::BicgstabWorkspace, A, b; kwargs...)
workspace = bicgstab!(workspace::BicgstabWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of bicgstab.
See BicgstabWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :bicgstab to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
DIOM
Krylov.diom — Function(x, stats) = diom(A, b::AbstractVector{FC};
memory::Int=20, M=I, N=I, ldiv::Bool=false,
reorthogonalization::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = diom(A, b, x0::AbstractVector; kwargs...)DIOM can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the consistent linear system Ax = b of size n 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.
An advantage of DIOM is that non-Hermitian or Hermitian indefinite or both non-Hermitian and indefinite systems of linear equations can be handled by this single algorithm.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :diom.
For an in-place variant that reuses memory across solves, see diom!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
memory: the number of most recent vectors of the Krylov basis against which to orthogonalize a new vector;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;reorthogonalization: reorthogonalize the new vectors of the Krylov basis against thememorymost recent vectors;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = diom!(workspace::DiomWorkspace, A, b; kwargs...)
workspace = diom!(workspace::DiomWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of diom. The keyword argument memory is the only exception. It is only supported by diom and is required to create a DiomWorkspace. It cannot be changed later.
See DiomWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :diom to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
FOM
Krylov.fom — Function(x, stats) = fom(A, b::AbstractVector{FC};
memory::Int=20, M=I, N=I, ldiv::Bool=false,
restart::Bool=false, reorthogonalization::Bool=false,
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = fom(A, b, x0::AbstractVector; kwargs...)FOM can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the linear system Ax = b of size n using FOM.
FOM algorithm is based on the Arnoldi process and a Galerkin condition.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :fom.
For an in-place variant that reuses memory across solves, see fom!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
memory: ifrestart = true, the restarted version FOM(k) is used withk = memory. Ifrestart = false, the parametermemoryshould be used as a hint of the number of iterations to limit dynamic memory allocations. Additional storage will be allocated if the number of iterations exceedsmemory;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;restart: restart the method aftermemoryiterations;reorthogonalization: reorthogonalize the new vectors of the Krylov basis against all previous vectors;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
Reference
- Y. Saad, Krylov subspace methods for solving unsymmetric linear systems, Mathematics of computation, Vol. 37(155), pp. 105–126, 1981.
Krylov.fom! — Functionworkspace = fom!(workspace::FomWorkspace, A, b; kwargs...)
workspace = fom!(workspace::FomWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of fom. The keyword argument memory is the only exception. It is only supported by fom and is required to create a FomWorkspace. It cannot be changed later.
See FomWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :fom to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
DQGMRES
Krylov.dqgmres — Function(x, stats) = dqgmres(A, b::AbstractVector{FC};
memory::Int=20, M=I, N=I, ldiv::Bool=false,
reorthogonalization::Bool=false, atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = dqgmres(A, b, x0::AbstractVector; kwargs...)DQGMRES can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the consistent linear system Ax = b of size n 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :dqgmres.
For an in-place variant that reuses memory across solves, see dqgmres!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
memory: the number of most recent vectors of the Krylov basis against which to orthogonalize a new vector;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;reorthogonalization: reorthogonalize the new vectors of the Krylov basis against thememorymost recent vectors;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = dqgmres!(workspace::DqgmresWorkspace, A, b; kwargs...)
workspace = dqgmres!(workspace::DqgmresWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of dqgmres. The keyword argument memory is the only exception. It is only supported by dqgmres and is required to create a DqgmresWorkspace. It cannot be changed later.
See DqgmresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :dqgmres to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
GMRES
Krylov.gmres — Function(x, stats) = gmres(A, b::AbstractVector{FC};
memory::Int=20, M=I, N=I, ldiv::Bool=false,
restart::Bool=false, reorthogonalization::Bool=false,
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = gmres(A, b, x0::AbstractVector; kwargs...)GMRES can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the linear system Ax = b of size n using GMRES.
GMRES algorithm is based on the Arnoldi process and computes a sequence of approximate solutions with the minimum residual.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :gmres.
For an in-place variant that reuses memory across solves, see gmres!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
memory: ifrestart = true, the restarted version GMRES(k) is used withk = memory. Ifrestart = false, the parametermemoryshould be used as a hint of the number of iterations to limit dynamic memory allocations. Additional storage will be allocated if the number of iterations exceedsmemory;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;restart: restart the method aftermemoryiterations;reorthogonalization: reorthogonalize the new vectors of the Krylov basis against all previous vectors;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Functionworkspace = gmres!(workspace::GmresWorkspace, A, b; kwargs...)
workspace = gmres!(workspace::GmresWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of gmres. The keyword argument memory is the only exception. It is only supported by gmres and is required to create a GmresWorkspace. It cannot be changed later.
See GmresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :gmres to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
FGMRES
Krylov.fgmres — Function(x, stats) = fgmres(A, b::AbstractVector{FC};
memory::Int=20, M=I, N=I, ldiv::Bool=false,
restart::Bool=false, reorthogonalization::Bool=false,
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(x, stats) = fgmres(A, b, x0::AbstractVector; kwargs...)FGMRES can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
Solve the linear system Ax = b of size n 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :fgmres.
For an in-place variant that reuses memory across solves, see fgmres!.
Input arguments
A: a linear operator that models a matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
memory: ifrestart = true, the restarted version FGMRES(k) is used withk = memory. Ifrestart = false, the parametermemoryshould be used as a hint of the number of iterations to limit dynamic memory allocations. Additional storage will be allocated if the number of iterations exceedsmemory;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;restart: restart the method aftermemoryiterations;reorthogonalization: reorthogonalize the new vectors of the Krylov basis against all previous vectors;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2n;timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
x: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
Reference
- Y. Saad, A Flexible Inner-Outer Preconditioned GMRES Algorithm, SIAM Journal on Scientific Computing, Vol. 14(2), pp. 461–469, 1993.
Krylov.fgmres! — Functionworkspace = fgmres!(workspace::FgmresWorkspace, A, b; kwargs...)
workspace = fgmres!(workspace::FgmresWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of fgmres. The keyword argument memory is the only exception. It is only supported by fgmres and is required to create a FgmresWorkspace. It cannot be changed later.
See FgmresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :fgmres to allocate the workspace, and krylov_solve! to run the Krylov method in-place.