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=solver->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.
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 n that represents an initial guess of the solution x.
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!
;radius
: add the trust-region constraint ‖x‖ ≤radius
ifradius > 0
. Useful to compute a step in a trust-region method for optimization;linesearch
: iftrue
, indicate that the solution is to be used in an inexact Newton method with linesearch. If negative curvature is detected at iteration k > 0, the solution of iteration k-1 is returned. If negative curvature is detected at iteration 0, the right-hand side is returned (i.e., the negative gradient);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 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 length n;stats
: statistics collected on the run in aSimpleStats
structure.
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!
— Functionsolver = cg!(solver::CgSolver, A, b; kwargs...)
solver = cg!(solver::CgSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of cg
.
See CgSolver
for more details about the solver
.
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=solver->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.
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 n that represents an initial guess of the solution x.
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!
;radius
: add the trust-region constraint ‖x‖ ≤radius
ifradius > 0
. Useful to compute a step in a trust-region method for optimization;linesearch
: iftrue
, indicate that the solution is to be used in an inexact Newton method with linesearch. If negative curvature is detected at iteration k > 0, the solution of iteration k-1 is returned. If negative curvature is detected at iteration 0, the right-hand side is returned (i.e., the negative gradient);γ
: 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 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 length n;stats
: statistics collected on the run in aSimpleStats
structure.
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!
— Functionsolver = cr!(solver::CrSolver, A, b; kwargs...)
solver = cr!(solver::CrSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of cr
.
See CrSolver
for more details about the solver
.
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=solver->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.
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 n that represents an initial guess of the solution x.
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!
;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 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 length n;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.car!
— Functionsolver = car!(solver::CarSolver, A, b; kwargs...)
solver = car!(solver::CarSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of car
.
See CarSolver
for more details about the solver
.
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=solver->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.
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 n that represents an initial guess of the solution x.
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!
;check_curvature
: iftrue
, check that the curvature of the quadratic along the search direction is positive, and abort if not, unlesslinesearch
is alsotrue
;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 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 length n;stats
: statistics collected on the run in aLanczosStats
structure.
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!
— Functionsolver = cg_lanczos!(solver::CgLanczosSolver, A, b; kwargs...)
solver = cg_lanczos!(solver::CgLanczosSolver, A, b, x0; kwargs...)
where kwargs
are keyword arguments of cg_lanczos
.
See CgLanczosSolver
for more details about the solver
.
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=solver->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.
Input arguments
A
: a linear operator that models a Hermitian matrix of dimension n;b
: a vector of length n;shifts
: a vector of length p.
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!
;check_curvature
: iftrue
, check that the curvature of the quadratic along the search direction is positive, and abort if not, unlesslinesearch
is alsotrue
;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 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 vector of p dense vectors, each one of length n;stats
: statistics collected on the run in aLanczosShiftStats
structure.
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!
— Functionsolver = cg_lanczos_shift!(solver::CgLanczosShiftSolver, A, b, shifts; kwargs...)
where kwargs
are keyword arguments of cg_lanczos_shift
.
See CgLanczosShiftSolver
for more details about the solver
.