SYMMLQ
Krylov.symmlq
— Function(x, stats) = symmlq(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false, window::Int=5,
transfer_to_cg::Bool=true, λ::T=zero(T),
λest::T=zero(T), etol::T=√eps(T),
conlim::T=1/√eps(T), atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=solver->false, iostream::IO=kstdout)
T
is an AbstractFloat
such as Float32
, Float64
or BigFloat
. FC
is T
or Complex{T}
.
(x, stats) = symmlq(A, b, x0::AbstractVector; kwargs...)
SYMMLQ can be warm-started from an initial guess x0
where kwargs
are the same keyword arguments as above
Solve the shifted linear system
(A + λI) x = b
of size n using the SYMMLQ method, where λ is a shift parameter, and A is Hermitian.
SYMMLQ produces monotonic errors ‖x* - x‖₂.
Input arguments
A
: a linear operator that models a Hermitian matrix of dimensionn
;b
: a vector of lengthn
.
Optional argument
x0
: a vector of lengthn
that represents an initial guess of the solutionx
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning;ldiv
: define whether the preconditioner usesldiv!
ormul!
;window
: number of iterations used to accumulate a lower bound on the error;transfer_to_cg
: transfer from the SYMMLQ point to the CG point, when it exists. The transfer is based on the residual norm;λ
: regularization parameter;λest
: positive strict lower bound on the smallest eigenvalueλₘᵢₙ
when solving a positive-definite system, such asλest = (1-10⁻⁷)λₘᵢₙ
;atol
: absolute stopping tolerance based on the residual norm;rtol
: relative stopping tolerance based on the residual norm;etol
: stopping tolerance based on the lower bound on the error;conlim
: limit on the estimated condition number ofA
beyond which the solution will be abandoned;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 everyverbose
iterations;history
: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback
: function or functor called ascallback(solver)
that returnstrue
if the Krylov method should terminate, andfalse
otherwise;iostream
: stream to which output is logged.
Output arguments
x
: a dense vector of lengthn
;stats
: statistics collected on the run in aSymmlqStats
structure.
Reference
- C. C. Paige and M. A. Saunders, Solution of Sparse Indefinite Systems of Linear Equations, SIAM Journal on Numerical Analysis, 12(4), pp. 617–629, 1975.
Krylov.symmlq!
— Functionsolver = symmlq!(solver::SymmlqSolver, A, b; kwargs...)
solver = symmlq!(solver::SymmlqSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of symmlq
.
See SymmlqSolver
for more details about the solver
.
MINRES
Krylov.minres
— Function(x, stats) = minres(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false, window::Int=5,
λ::T=zero(T), atol::T=√eps(T),
rtol::T=√eps(T), etol::T=√eps(T),
conlim::T=1/√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=solver->false, iostream::IO=kstdout)
T
is an AbstractFloat
such as Float32
, Float64
or BigFloat
. FC
is T
or Complex{T}
.
(x, stats) = minres(A, b, x0::AbstractVector; kwargs...)
MINRES can be warm-started from an initial guess x0
where kwargs
are the same keyword arguments as above.
Solve the shifted linear least-squares problem
minimize ‖b - (A + λI)x‖₂²
or the shifted linear system
(A + λI) x = b
of size n using the MINRES method, where λ ≥ 0 is a shift parameter, where A is Hermitian.
MINRES is formally equivalent to applying CR to Ax=b when A is positive definite, but is typically more stable and also applies to the case where A is indefinite.
MINRES produces monotonic residuals ‖r‖₂ and optimality residuals ‖Aᴴr‖₂.
Input arguments
A
: a linear operator that models a Hermitian matrix of dimensionn
;b
: a vector of lengthn
.
Optional argument
x0
: a vector of lengthn
that represents an initial guess of the solutionx
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning;ldiv
: define whether the preconditioner usesldiv!
ormul!
;window
: number of iterations used to accumulate a lower bound on the error;λ
: regularization parameter;atol
: absolute stopping tolerance based on the residual norm;rtol
: relative stopping tolerance based on the residual norm;etol
: stopping tolerance based on the lower bound on the error;conlim
: limit on the estimated condition number ofA
beyond which the solution will be abandoned;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 everyverbose
iterations;history
: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback
: function or functor called ascallback(solver)
that returnstrue
if the Krylov method should terminate, andfalse
otherwise;iostream
: stream to which output is logged.
Output arguments
x
: a dense vector of lengthn
;stats
: statistics collected on the run in aSimpleStats
structure.
Reference
- C. C. Paige and M. A. Saunders, Solution of Sparse Indefinite Systems of Linear Equations, SIAM Journal on Numerical Analysis, 12(4), pp. 617–629, 1975.
Krylov.minres!
— Functionsolver = minres!(solver::MinresSolver, A, b; kwargs...)
solver = minres!(solver::MinresSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of minres
.
See MinresSolver
for more details about the solver
.
MINRES-QLP
Krylov.minres_qlp
— Function(x, stats) = minres_qlp(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false, Artol::T=√eps(T),
λ::T=zero(T), atol::T=√eps(T),
rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=solver->false, iostream::IO=kstdout)
T
is an AbstractFloat
such as Float32
, Float64
or BigFloat
. FC
is T
or Complex{T}
.
(x, stats) = minres_qlp(A, b, x0::AbstractVector; kwargs...)
MINRES-QLP can be warm-started from an initial guess x0
where kwargs
are the same keyword arguments as above.
MINRES-QLP is the only method based on the Lanczos process that returns the minimum-norm solution on singular inconsistent systems (A + λI)x = b of size n, where λ is a shift parameter. It is significantly more complex but can be more reliable than MINRES when A is ill-conditioned.
M also indicates the weighted norm in which residuals are measured.
Input arguments
A
: a linear operator that models a Hermitian matrix of dimensionn
;b
: a vector of lengthn
.
Optional argument
x0
: a vector of lengthn
that represents an initial guess of the solutionx
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning;ldiv
: define whether the preconditioner usesldiv!
ormul!
;λ
: regularization parameter;atol
: absolute stopping tolerance based on the residual norm;rtol
: relative stopping tolerance based on the residual norm;Artol
: relative stopping tolerance based on the Aᴴ-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 everyverbose
iterations;history
: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback
: function or functor called ascallback(solver)
that returnstrue
if the Krylov method should terminate, andfalse
otherwise;iostream
: stream to which output is logged.
Output arguments
x
: a dense vector of lengthn
;stats
: statistics collected on the run in aSimpleStats
structure.
References
- S.-C. T. Choi, Iterative methods for singular linear equations and least-squares problems, Ph.D. thesis, ICME, Stanford University, 2006.
- S.-C. T. Choi, C. C. Paige and M. A. Saunders, MINRES-QLP: A Krylov subspace method for indefinite or singular symmetric systems, SIAM Journal on Scientific Computing, Vol. 33(4), pp. 1810–1836, 2011.
- S.-C. T. Choi and M. A. Saunders, Algorithm 937: MINRES-QLP for symmetric and Hermitian linear equations and least-squares problems, ACM Transactions on Mathematical Software, 40(2), pp. 1–12, 2014.
Krylov.minres_qlp!
— Functionsolver = minres_qlp!(solver::MinresQlpSolver, A, b; kwargs...)
solver = minres_qlp!(solver::MinresQlpSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of minres_qlp
.
See MinresQlpSolver
for more details about the solver
.
MINARES
Krylov.minares
— Function(x, stats) = minares(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false,
λ::T = zero(T), atol::T=√eps(T),
rtol::T=√eps(T), Artol::T = √eps(T),
itmax::Int=0, timemax::Float64=Inf,
verbose::Int=0, history::Bool=false,
callback=solver->false, iostream::IO=kstdout)
T
is an AbstractFloat
such as Float32
, Float64
or BigFloat
. FC
is T
or Complex{T}
.
(x, stats) = minares(A, b, x0::AbstractVector; kwargs...)
MINARES can be warm-started from an initial guess x0
where kwargs
are the same keyword arguments as above.
MINARES solves the Hermitian linear system Ax = b of size n. MINARES minimizes ‖Arₖ‖₂ when M = Iₙ and ‖AMrₖ‖M otherwise. The estimates computed every iteration are ‖Mrₖ‖₂ and ‖AMrₖ‖M.
Input arguments
A
: a linear operator that models a Hermitian positive definite matrix of dimensionn
;b
: a vector of lengthn
.
Optional argument
x0
: a vector of lengthn
that represents an initial guess of the solutionx
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning;ldiv
: define whether the preconditioner usesldiv!
ormul!
;λ
: regularization parameter;atol
: absolute stopping tolerance based on the residual norm;rtol
: relative stopping tolerance based on the residual norm;Artol
: relative stopping tolerance based on the Aᴴ-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 everyverbose
iterations;history
: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;callback
: function or functor called ascallback(solver)
that returnstrue
if the Krylov method should terminate, andfalse
otherwise;iostream
: stream to which output is logged.
Output arguments
x
: a dense vector of lengthn
;stats
: statistics collected on the run in aSimpleStats
structure.
Reference
- A. Montoison, D. Orban and M. A. Saunders, MinAres: An Iterative Solver for Symmetric Linear Systems, Cahier du GERAD G-2023-40, GERAD, Montréal, 2023.
Krylov.minares!
— Functionsolver = minares!(solver::MinaresSolver, A, b; kwargs...)
solver = minares!(solver::MinaresSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of minares
.
See MinaresSolver
for more details about the solver
.