CG
Krylov.cg — Function(x, stats) = cg(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false, radius::T=zero(T),
linesearch::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) = cg(A, b, x0::AbstractVector; kwargs...)CG can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
The conjugate gradient method to solve the Hermitian linear system Ax = b of size n.
The method does not abort if A is not definite. 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 = :cg.
For an in-place variant that reuses memory across solves, see cg!.
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 lengthnthat represents an initial guess of the solutionx.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioner usesldiv!ormul!;radius: add the trust-region constraint ‖x‖ ≤radiusifradius > 0. Useful to compute a step in a trust-region method for optimization.- If 'radius' > 0, and nonpositive curvature is detected along the current search direction, we take the step to the trust-region boundary, the search direction is stored in
stats.npc_dir, andstats.npcCountis set to 1.
- If 'radius' > 0, and nonpositive curvature is detected along the current search direction, we take the step to the trust-region boundary, the search direction is stored in
linesearch: whentrue, assume that CG is used within an inexact Newton method with line search.- If nonpositive curvature occurs at k = 0, the solver instead 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, andstats.npcCountis set to 1.; - If nonpositive curvature is detected at iteration k > 0, the method rolls back to the solution from iteration k – 1. The search direction computed at iteration k is stored in
stats.npc_dir, andstats.npcCountis set to 1.
- If nonpositive curvature occurs at k = 0, the solver instead takes the right-hand side (i.e., the preconditioned negative gradient) as the current solution. The same search direction is returned in
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
- M. R. Hestenes and E. Stiefel, Methods of conjugate gradients for solving linear systems, Journal of Research of the National Bureau of Standards, 49(6), pp. 409–436, 1952.
Krylov.cg! — Functionworkspace = cg!(workspace::CgWorkspace, A, b; kwargs...)
workspace = cg!(workspace::CgWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of cg.
See CgWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :cg to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
CR
Krylov.cr — Function(x, stats) = cr(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false, radius::T=zero(T),
linesearch::Bool=false, γ::T=√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) = cr(A, b, x0::AbstractVector; kwargs...)CR can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
A truncated version of Stiefel’s Conjugate Residual method to solve the Hermitian linear system Ax = b of size n or the least-squares problem min ‖b - Ax‖ if A is singular. The matrix A must be Hermitian semi-definite. 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 = :cr.
For an in-place variant that reuses memory across solves, see cr!.
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 lengthnthat represents an initial guess of the solutionx.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioner usesldiv!ormul!;radius: add the trust-region constraint ‖x‖ ≤radiusifradius > 0. Useful to compute a step in a trust-region method for optimization; Ifradius > 0and nonpositive curvature is detected, the behavior depends on the iteration and follow similar logic as linesearch;linesearch: iftrueand nonpositive curvature is detected, the behavior depends on the iteration:
– at iteration k = 0, return the preconditioned initial search direction in workspace.npc_dir; – at iteration k > 0,
- if the residual from iteration k-1 is a nonpositive curvature direction but
workspace.p, the search direction at iteration k, is not, the residual is stored instats.npc_dirandstats.npcCountis set to 1; - if
workspace.pis a nonpositive curvature direction but the residual is not,workspace.pis copied intostats.npc_dirandstats.npcCountis set to 1; - if both are nonpositive curvature directions, the residual is stored in
stats.npc_dirandstats.npcCountis set to 2. γ: tolerance to determine that the curvature of the quadratic model is nonpositive;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
- M. R. Hestenes and E. Stiefel, Methods of conjugate gradients for solving linear systems, Journal of Research of the National Bureau of Standards, 49(6), pp. 409–436, 1952.
- E. Stiefel, Relaxationsmethoden bester Strategie zur Losung linearer Gleichungssysteme, Commentarii Mathematici Helvetici, 29(1), pp. 157–179, 1955.
- D. G. Luenberger, The conjugate residual method for constrained minimization problems, SIAM Journal on Numerical Analysis, 7(3), pp. 390–398, 1970.
- 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.cr! — Functionworkspace = cr!(workspace::CrWorkspace, A, b; kwargs...)
workspace = cr!(workspace::CrWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of cr.
See CrWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :cr to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
CAR
Krylov.car — Function(x, stats) = car(A, b::AbstractVector{FC};
M=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) = car(A, b, x0::AbstractVector; kwargs...)CAR can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
CAR solves the Hermitian and positive definite linear system Ax = b of size n. CAR 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 = :car.
For an in-place variant that reuses memory across solves, see car!.
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 lengthnthat represents an initial guess of the solutionx.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioner usesldiv!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
- 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.car! — Functionworkspace = car!(workspace::CarWorkspace, A, b; kwargs...)
workspace = car!(workspace::CarWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of car.
See CarWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :car to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
CG-LANCZOS
Krylov.cg_lanczos — Function(x, stats) = cg_lanczos(A, b::AbstractVector{FC};
M=I, ldiv::Bool=false,
check_curvature::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) = cg_lanczos(A, b, x0::AbstractVector; kwargs...)CG-LANCZOS can be warm-started from an initial guess x0 where kwargs are the same keyword arguments as above.
The Lanczos version of the conjugate gradient method to solve the Hermitian linear system Ax = b of size n.
The method does not abort if A is not definite.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :cg_lanczos.
For an in-place variant that reuses memory across solves, see cg_lanczos!.
Input arguments
A: a linear operator that models a Hermitian matrix of dimensionn;b: a vector of lengthn.
Optional argument
x0: a vector of lengthnthat represents an initial guess of the solutionx.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioner usesldiv!ormul!;check_curvature: iftrue, check that the curvature of the quadratic along the search direction is positive, and abort if not;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 aLanczosStatsstructure.
References
- A. Frommer and P. Maass, Fast CG-Based Methods for Tikhonov-Phillips Regularization, SIAM Journal on Scientific Computing, 20(5), pp. 1831–1850, 1999.
- 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.cg_lanczos! — Functionworkspace = cg_lanczos!(workspace::CgLanczosWorkspace, A, b; kwargs...)
workspace = cg_lanczos!(workspace::CgLanczosWorkspace, A, b, x0; kwargs...)In these calls, kwargs are keyword arguments of cg_lanczos.
See CgLanczosWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :cg_lanczos to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
CG-LANCZOS-SHIFT
Krylov.cg_lanczos_shift — Function(x, stats) = cg_lanczos_shift(A, b::AbstractVector{FC}, shifts::AbstractVector{T};
M=I, ldiv::Bool=false,
check_curvature::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}.
The Lanczos version of the conjugate gradient method to solve a family of shifted systems
(A + αI) x = b (α = α₁, ..., αₚ)of size n. The method does not abort if A + αI is not definite.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :cg_lanczos_shift.
For an in-place variant that reuses memory across solves, see cg_lanczos_shift!.
Input arguments
A: a linear operator that models a Hermitian matrix of dimensionn;b: a vector of lengthn;shifts: a vector of lengthp.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioner usesldiv!ormul!;check_curvature: iftrue, check that the curvature of the quadratic along the search direction is positive, and abort if not;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 vector ofpdense vectors, each one of lengthn;stats: statistics collected on the run in aLanczosShiftStatsstructure.
References
- A. Frommer and P. Maass, Fast CG-Based Methods for Tikhonov-Phillips Regularization, SIAM Journal on Scientific Computing, 20(5), pp. 1831–1850, 1999.
- 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.cg_lanczos_shift! — Functionworkspace = cg_lanczos_shift!(workspace::CgLanczosShiftWorkspace, A, b, shifts; kwargs...)In this call, kwargs are keyword arguments of cg_lanczos_shift.
See CgLanczosShiftWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :cg_lanczos_shift to allocate the workspace, and krylov_solve! to run the Krylov method in-place.