TriCG
Krylov.tricg
— Function(x, y, stats) = tricg(A, b::AbstractVector{FC}, c::AbstractVector{FC};
M=I, N=I, ldiv::Bool=false,
spd::Bool=false, snd::Bool=false,
flip::Bool=false, τ::T=one(T),
ν::T=-one(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, y, stats) = tricg(A, b, c, x0::AbstractVector, y0::AbstractVector; kwargs...)
TriCG can be warm-started from initial guesses x0
and y0
where kwargs
are the same keyword arguments as above.
Given a matrix A
of dimension m × n, TriCG solves the Hermitian linear system
[ τE A ] [ x ] = [ b ]
[ Aᴴ νF ] [ y ] [ c ],
of size (n+m) × (n+m) 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 Hermitian and quasi-definite linear systems with τ = 1 and ν = -1.
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.
Input arguments
A
: a linear operator that models a matrix of dimensionm × n
;b
: a vector of lengthm
;c
: a vector of lengthn
.
Optional arguments
x0
: a vector of lengthm
that represents an initial guess of the solutionx
;y0
: a vector of lengthn
that represents an initial guess of the solutiony
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizem
used for centered preconditioning of the partitioned system;N
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning of the partitioned system;ldiv
: define whether the preconditioners useldiv!
ormul!
;spd
: iftrue
, setτ = 1
andν = 1
for Hermitian and positive-definite linear system;snd
: iftrue
, setτ = -1
andν = -1
for Hermitian and negative-definite linear systems;flip
: iftrue
, setτ = -1
andν = 1
for another known variant of Hermitian quasi-definite systems;τ
andν
: diagonal scaling factors of the partitioned Hermitian linear system;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 tom+n
;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 lengthm
;y
: a dense vector of lengthn
;stats
: statistics collected on the run in aSimpleStats
structure.
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, ldiv::Bool=false,
spd::Bool=false, snd::Bool=false,
flip::Bool=false, sp::Bool=false,
τ::T=one(T), ν::T=-one(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, y, stats) = trimr(A, b, c, x0::AbstractVector, y0::AbstractVector; kwargs...)
TriMR can be warm-started from initial guesses x0
and y0
where kwargs
are the same keyword arguments as above.
Given a matrix A
of dimension m × n, TriMR solves the symmetric linear system
[ τE A ] [ x ] = [ b ]
[ Aᴴ νF ] [ y ] [ c ],
of size (n+m) × (n+m) 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.
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.
Input arguments
A
: a linear operator that models a matrix of dimensionm × n
;b
: a vector of lengthm
;c
: a vector of lengthn
.
Optional arguments
x0
: a vector of lengthm
that represents an initial guess of the solutionx
;y0
: a vector of lengthn
that represents an initial guess of the solutiony
.
Keyword arguments
M
: linear operator that models a Hermitian positive-definite matrix of sizem
used for centered preconditioning of the partitioned system;N
: linear operator that models a Hermitian positive-definite matrix of sizen
used for centered preconditioning of the partitioned system;ldiv
: define whether the preconditioners useldiv!
ormul!
;spd
: iftrue
, setτ = 1
andν = 1
for Hermitian and positive-definite linear system;snd
: iftrue
, setτ = -1
andν = -1
for Hermitian and negative-definite linear systems;flip
: iftrue
, setτ = -1
andν = 1
for another known variant of Hermitian quasi-definite systems;sp
: iftrue
, setτ = 1
andν = 0
for saddle-point systems;τ
andν
: diagonal scaling factors of the partitioned Hermitian linear system;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 tom+n
;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 lengthm
;y
: a dense vector of lengthn
;stats
: statistics collected on the run in aSimpleStats
structure.
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
.