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=workspace->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 = bof size n using the SYMMLQ method, where λ is a shift parameter, and A is Hermitian.
SYMMLQ produces monotonic errors ‖x* - x‖₂.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :symmlq.
For an in-place variant that reuses memory across solves, see symmlq!.
Input arguments
- A: a linear operator that models a Hermitian matrix of dimension- n;
- b: a vector of length- n.
Optional argument
- x0: a vector of length- nthat represents an initial guess of the solution- x.
Keyword arguments
- M: linear operator that models a Hermitian positive-definite matrix of size- nused for centered preconditioning;
- ldiv: define whether the preconditioner uses- ldiv!or- mul!;
- 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 of- Abeyond which the solution will be abandoned;
- itmax: the maximum number of iterations. If- itmax=0, the default number of iterations is set to- 2n;
- timemax: the time limit in seconds;
- verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every- verboseiterations;
- history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;
- callback: function or functor called as- callback(workspace)that returns- trueif the Krylov method should terminate, and- falseotherwise;
- iostream: stream to which output is logged.
Output arguments
- x: a dense vector of length- n;
- stats: statistics collected on the run in a- SymmlqStatsstructure.
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! — Functionworkspace = symmlq!(workspace::SymmlqWorkspace, A, b; kwargs...)
workspace = symmlq!(workspace::SymmlqWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of symmlq.
See SymmlqWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :symmlq to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
MINRES
Krylov.minres — Function(x, stats) = minres(A, b::AbstractVector{FC};
                    M=I, ldiv::Bool=false, window::Int=5,
                    linesearch::Bool=false, λ::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=workspace->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 = bof 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‖₂.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :minres.
For an in-place variant that reuses memory across solves, see minres!.
Input arguments
- A: a linear operator that models a Hermitian matrix of dimension- n;
- b: a vector of length- n.
Optional argument
- x0: a vector of length- nthat represents an initial guess of the solution- x.
Keyword arguments
- M: linear operator that models a Hermitian positive-definite matrix of size- nused for centered preconditioning;
- ldiv: define whether the preconditioner uses- ldiv!or- mul!;
- window: number of iterations used to accumulate a lower bound on the error;
- linesearch: if- true, indicate that the solution is to be used in an inexact Newton method with linesearch. If- trueand nonpositive curvature is detected, the behavior depends on the iteration:
– at iteration k = 1, the solver takes the right-hand side (i.e., the preconditioned negative gradient) as the current solution. The same search direction is returned in workspace.npc_dir, and stats.npcCount is set to 1;  – at iteration k > 1, the solver returns the solution from iteration k – 1,
- if the residual from iteration k is a nonpositive curvature direction but the search direction at iteration k, is not, the residual is stored in stats.npc_dirandstats.npcCountis set to 1;
- if both are nonpositive curvature directions, the residual is stored in stats.npc_dir, the search direction is stored inworkspace.w1, andstats.npcCountis set to 2. (Note that the MINRES solver starts at iteration 1, so the first iteration is k = 1);
- λ: 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 of- Abeyond which the solution will be abandoned;
- itmax: the maximum number of iterations. If- itmax=0, the default number of iterations is set to- 2n;
- timemax: the time limit in seconds;
- verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every- verboseiterations;
- history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;
- callback: function or functor called as- callback(workspace)that returns- trueif the Krylov method should terminate, and- falseotherwise;
- iostream: stream to which output is logged.
Output arguments
- x: a dense vector of length- n;
- stats: statistics collected on the run in a- SimpleStatsstructure.
References
- 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. 
- Y. Liu and F. Roosta, MINRES: From Negative Curvature Detection to Monotonicity Properties, SIAM Journal on Optimization, 32(4), pp. 2636–2661, 2022. 
- M-A. Dahito and D. Orban, The Conjugate Residual Method in Linesearch and Trust-Region Methods, SIAM Journal on Optimization, 29(3), pp. 1988–2025, 2019. 
Krylov.minres! — Functionworkspace = minres!(workspace::MinresWorkspace, A, b; kwargs...)
workspace = minres!(workspace::MinresWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of minres.
See MinresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :minres to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
MINRES-QLP
Krylov.minres_qlp — Function(x, stats) = minres_qlp(A, b::AbstractVector{FC};
                        M=I, ldiv::Bool=false, Artol::T=√eps(T),
                        linesearch::Bool=false, λ::T=zero(T), 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) = 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :minres_qlp.
For an in-place variant that reuses memory across solves, see minres_qlp!.
Input arguments
- A: a linear operator that models a Hermitian matrix of dimension- n;
- b: a vector of length- n.
Optional argument
- x0: a vector of length- nthat represents an initial guess of the solution- x.
Keyword arguments
- M: linear operator that models a Hermitian positive-definite matrix of size- nused for centered preconditioning;
- ldiv: define whether the preconditioner uses- ldiv!or- mul!;
- λ: 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;
- linesearch: if- true, indicate that the solution is to be used in an inexact Newton method with linesearch. If- trueand nonpositive curvature is detected, the behavior depends on the iteration:
– at iteration k = 1, the solver takes the right-hand side (i.e., the preconditioned negative gradient) as the current solution. The same search direction is returned in workspace.npc_dir, and stats.npcCount is set to 1;  – at iteration k > 1, the solver returns the solution from iteration k – 1, the residual from iteration k is a nonpositive curvature direction stored in stats.npc_dir and stats.npcCount is set to 1;
- itmax: the maximum number of iterations. If- itmax=0, the default number of iterations is set to- 2n;
- timemax: the time limit in seconds;
- verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every- verboseiterations;
- history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;
- callback: function or functor called as- callback(workspace)that returns- trueif the Krylov method should terminate, and- falseotherwise;
- iostream: stream to which output is logged.
Output arguments
- x: a dense vector of length- n;
- stats: statistics collected on the run in a- SimpleStatsstructure.
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.
- Y. Liu and F. Roosta, MINRES: From Negative Curvature Detection to Monotonicity Properties, SIAM Journal on Optimization, 32(4), pp. 2636–2661, 2022.
Krylov.minres_qlp! — Functionworkspace = minres_qlp!(workspace::MinresQlpWorkspace, A, b; kwargs...)
workspace = minres_qlp!(workspace::MinresQlpWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of minres_qlp.
See MinresQlpWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :minres_qlp to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
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=workspace->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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :minares.
For an in-place variant that reuses memory across solves, see minares!.
Input arguments
- A: a linear operator that models a Hermitian positive definite matrix of dimension- n;
- b: a vector of length- n.
Optional argument
- x0: a vector of length- nthat represents an initial guess of the solution- x.
Keyword arguments
- M: linear operator that models a Hermitian positive-definite matrix of size- nused for centered preconditioning;
- ldiv: define whether the preconditioner uses- ldiv!or- mul!;
- λ: 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. If- itmax=0, the default number of iterations is set to- 2n;
- timemax: the time limit in seconds;
- verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every- verboseiterations;
- history: collect additional statistics on the run such as residual norms, or Aᴴ-residual norms;
- callback: function or functor called as- callback(workspace)that returns- trueif the Krylov method should terminate, and- falseotherwise;
- iostream: stream to which output is logged.
Output arguments
- x: a dense vector of length- n;
- stats: statistics collected on the run in a- SimpleStatsstructure.
Reference
- A. Montoison, D. Orban and M. A. Saunders, MinAres: An Iterative Solver for Symmetric Linear Systems, SIAM Journal on Matrix Analysis and Applications, 46(1), pp. 509–529, 2025.
Krylov.minares! — Functionworkspace = minares!(workspace::MinaresWorkspace, A, b; kwargs...)
workspace = minares!(workspace::MinaresWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of minares.
See MinaresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :minares to allocate the workspace, and krylov_solve! to run the Krylov method in-place.