CG
Krylov.cg — Function(x, stats) = cg(A, b::AbstractVector{FC};
M=I, atol::T=√eps(T), rtol::T=√eps(T),
itmax::Int=0, radius::T=zero(T), linesearch::Bool=false,
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}.
The conjugate gradient method to solve the symmetric linear system Ax=b.
The method does not abort if A is not definite.
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.
If itmax=0, the default number of iterations is set to 2 * n, with n = length(b).
CG can be warm-started from an initial guess x0 with the method
(x, stats) = cg(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
- 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, atol::T=√eps(T), rtol::T=√eps(T), γ::T=√eps(T), itmax::Int=0,
radius::T=zero(T), verbose::Int=0, linesearch::Bool=false, 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}.
A truncated version of Stiefel’s Conjugate Residual method to solve the symmetric linear system Ax = b or the least-squares problem min ‖b - Ax‖. The matrix A must be positive semi-definite.
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.
In a linesearch context, 'linesearch' must be set to 'true'.
If itmax=0, the default number of iterations is set to 2 * n, with n = length(b).
CR can be warm-started from an initial guess x0 with the method
(x, stats) = cr(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
- 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.
- 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.
CG-LANCZOS
Krylov.cg_lanczos — Function(x, stats) = cg_lanczos(A, b::AbstractVector{FC};
M=I, atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
check_curvature::Bool=false, 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}.
The Lanczos version of the conjugate gradient method to solve the symmetric linear system
Ax = bThe method does not abort if A is not definite.
A preconditioner M may be provided in the form of a linear operator and is assumed to be hermitian and positive definite.
CG-LANCZOS can be warm-started from an initial guess x0 with the method
(x, stats) = cg_lanczos(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
- 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, atol::T=√eps(T), rtol::T=√eps(T),
itmax::Int=0, check_curvature::Bool=false,
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}.
The Lanczos version of the conjugate gradient method to solve a family of shifted systems
(A + αI) x = b (α = α₁, ..., αₙ)The method does not abort if A + αI is not definite.
A preconditioner M may be provided in the form of a linear operator and is assumed to be hermitian and positive definite.
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
Krylov.cg_lanczos_shift! — Functionsolver = cg_lanczos!(solver::CgLanczosShiftSolver, A, b, shifts; kwargs...)where kwargs are keyword arguments of cg_lanczos_shift.
See CgLanczosShiftSolver for more details about the solver.