TriCG
Krylov.tricg
— Function(x, y, stats) = tricg(A, b::AbstractVector{FC}, c::AbstractVector{FC};
M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
spd::Bool=false, snd::Bool=false, flip::Bool=false,
τ::T=one(T), ν::T=-one(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}
.
TriCG solves the symmetric linear system
[ τE A ] [ x ] = [ b ]
[ Aᵀ νF ] [ y ] [ c ],
where τ and ν are real numbers, E = M⁻¹ ≻ 0 and F = N⁻¹ ≻ 0. b
and c
must both be nonzero. TriCG could breakdown if τ = 0
or ν = 0
. It's recommended to use TriMR in these cases.
By default, TriCG solves symmetric and quasi-definite linear systems with τ = 1 and ν = -1. If flip = true
, TriCG solves another known variant of SQD systems where τ = -1 and ν = 1. If spd = true
, τ = ν = 1 and the associated symmetric and positive definite linear system is solved. If snd = true
, τ = ν = -1 and the associated symmetric and negative definite linear system is solved. τ
and ν
are also keyword arguments that can be directly modified for more specific problems.
TriCG is based on the preconditioned orthogonal tridiagonalization process and its relation with the preconditioned block-Lanczos process.
[ M 0 ]
[ 0 N ]
indicates the weighted norm in which residuals are measured. It's the Euclidean norm when M
and N
are identity operators.
TriCG stops when itmax
iterations are reached or when ‖rₖ‖ ≤ atol + ‖r₀‖ * rtol
. atol
is an absolute tolerance and rtol
is a relative tolerance.
Additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every verbose
iterations.
TriCG can be warm-started from initial guesses x0
and y0
with the method
(x, y, stats) = tricg(A, b, c, x0, y0; 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
- A. Montoison and D. Orban, TriCG and TriMR: Two Iterative Methods for Symmetric Quasi-Definite Systems, SIAM Journal on Scientific Computing, 43(4), pp. 2502–2525, 2021.
Krylov.tricg!
— Functionsolver = tricg!(solver::TricgSolver, A, b, c; kwargs...)
solver = tricg!(solver::TricgSolver, A, b, c, x0, y0; kwargs...)
where kwargs
are keyword arguments of tricg
.
See TricgSolver
for more details about the solver
.
TriMR
Krylov.trimr
— Function(x, y, stats) = trimr(A, b::AbstractVector{FC}, c::AbstractVector{FC};
M=I, N=I, atol::T=√eps(T), rtol::T=√eps(T),
spd::Bool=false, snd::Bool=false, flip::Bool=false, sp::Bool=false,
τ::T=one(T), ν::T=-one(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}
.
TriMR solves the symmetric linear system
[ τE A ] [ x ] = [ b ]
[ Aᵀ νF ] [ y ] [ c ],
where τ and ν are real numbers, E = M⁻¹ ≻ 0, F = N⁻¹ ≻ 0. b
and c
must both be nonzero. TriMR handles saddle-point systems (τ = 0
or ν = 0
) and adjoint systems (τ = 0
and ν = 0
) without any risk of breakdown.
By default, TriMR solves symmetric and quasi-definite linear systems with τ = 1 and ν = -1. If flip = true
, TriMR solves another known variant of SQD systems where τ = -1 and ν = 1. If spd = true
, τ = ν = 1 and the associated symmetric and positive definite linear system is solved. If snd = true
, τ = ν = -1 and the associated symmetric and negative definite linear system is solved. If sp = true
, τ = 1, ν = 0 and the associated saddle-point linear system is solved. τ
and ν
are also keyword arguments that can be directly modified for more specific problems.
TriMR is based on the preconditioned orthogonal tridiagonalization process and its relation with the preconditioned block-Lanczos process.
[ M 0 ]
[ 0 N ]
indicates the weighted norm in which residuals are measured. It's the Euclidean norm when M
and N
are identity operators.
TriMR stops when itmax
iterations are reached or when ‖rₖ‖ ≤ atol + ‖r₀‖ * rtol
. atol
is an absolute tolerance and rtol
is a relative tolerance.
Additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed every verbose
iterations.
TriMR can be warm-started from initial guesses x0
and y0
with the method
(x, y, stats) = trimr(A, b, c, x0, y0; 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
- A. Montoison and D. Orban, TriCG and TriMR: Two Iterative Methods for Symmetric Quasi-Definite Systems, SIAM Journal on Scientific Computing, 43(4), pp. 2502–2525, 2021.
Krylov.trimr!
— Functionsolver = trimr!(solver::TrimrSolver, A, b, c; kwargs...)
solver = trimr!(solver::TrimrSolver, A, b, c, x0, y0; kwargs...)
where kwargs
are keyword arguments of trimr
.
See TrimrSolver
for more details about the solver
.