SYMMLQ
Krylov.symmlq — Function(x, stats) = symmlq(A, b::AbstractVector{FC}; window::Int=0,
M=I, λ::T=zero(T), transfer_to_cg::Bool=true,
λest::T=zero(T), atol::T=√eps(T), rtol::T=√eps(T),
etol::T=√eps(T), itmax::Int=0, conlim::T=1/√eps(T),
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 shifted linear system
(A + λI) x = busing the SYMMLQ method, where λ is a shift parameter, and A is square and symmetric.
SYMMLQ produces monotonic errors ‖x*-x‖₂.
A preconditioner M may be provided in the form of a linear operator and is assumed to be symmetric and positive definite.
SYMMLQ can be warm-started from an initial guess x0 with the method
(x, stats) = symmlq(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
- 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, λ::T=zero(T), atol::T=√eps(T)/100,
rtol::T=√eps(T)/100, ratol :: T=zero(T),
rrtol :: T=zero(T), etol::T=√eps(T),
window::Int=5, itmax::Int=0,
conlim::T=1/√eps(T), 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 shifted linear least-squares problem
minimize ‖b - (A + λI)x‖₂²or the shifted linear system
(A + λI) x = busing the MINRES method, where λ ≥ 0 is a shift parameter, where A is square and symmetric.
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‖₂.
A preconditioner M may be provided in the form of a linear operator and is assumed to be symmetric and positive definite.
MINRES can be warm-started from an initial guess x0 with the method
(x, stats) = minres(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
- 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, atol::T=√eps(T), rtol::T=√eps(T),
ctol::T=√eps(T), λ::T=zero(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}.
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, where λ is a shift parameter. It is significantly more complex but can be more reliable than MINRES when A is ill-conditioned.
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.
MINRES-QLP can be warm-started from an initial guess x0 with the method
(x, stats) = minres_qlp(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
- 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.