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=workspace->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. 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.
The matrix
[ M 0 ]
[ 0 N ]indicates the weighted norm in which residuals are measured, here denoted ‖·‖. 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :tricg.
For an in-place variant that reuses memory across solves, see tricg!.
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 lengthmthat represents an initial guess of the solutionx;y0: a vector of lengthnthat represents an initial guess of the solutiony.
Warm-starting is supported only when M and N are either I or the corresponding coefficient (τ or ν) is zero.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizemused for centered preconditioning of the partitioned system;N: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning of the partitioned system;ldiv: define whether the preconditioners useldiv!ormul!;spd: iftrue, setτ = 1andν = 1for Hermitian and positive-definite linear system;snd: iftrue, setτ = -1andν = -1for Hermitian and negative-definite linear systems;flip: iftrue, setτ = -1andν = 1for 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 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 lengthm;y: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Function
workspace = tricg!(workspace::TricgWorkspace, A, b, c; kwargs...)
workspace = tricg!(workspace::TricgWorkspace, A, b, c, x0, y0; kwargs...)In these calls, kwargs are keyword arguments of tricg.
See TricgWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :tricg to allocate the workspace, and krylov_solve! to run the Krylov method in-place.
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=workspace->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. 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.
The matrix
[ M 0 ]
[ 0 N ]indicates the weighted norm in which residuals are measured, here denoted ‖·‖. 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.
Interface
To easily switch between Krylov methods, use the generic interface krylov_solve with method = :trimr.
For an in-place variant that reuses memory across solves, see trimr!.
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 lengthmthat represents an initial guess of the solutionx;y0: a vector of lengthnthat represents an initial guess of the solutiony.
Warm-starting is supported only when M and N are either I or the corresponding coefficient (τ or ν) is zero.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizemused for centered preconditioning of the partitioned system;N: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning of the partitioned system;ldiv: define whether the preconditioners useldiv!ormul!;spd: iftrue, setτ = 1andν = 1for Hermitian and positive-definite linear system;snd: iftrue, setτ = -1andν = -1for Hermitian and negative-definite linear systems;flip: iftrue, setτ = -1andν = 1for another known variant of Hermitian quasi-definite systems;sp: iftrue, setτ = 1andν = 0for 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 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 lengthm;y: a dense vector of lengthn;stats: statistics collected on the run in aSimpleStatsstructure.
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! — Function
workspace = trimr!(workspace::TrimrWorkspace, A, b, c; kwargs...)
workspace = trimr!(workspace::TrimrWorkspace, A, b, c, x0, y0; kwargs...)In these calls, kwargs are keyword arguments of trimr.
See TrimrWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :trimr to allocate the workspace, and krylov_solve! to run the Krylov method in-place.