MINRES
MINRES is a short-recurrence version of GMRES for solving $Ax = b$ approximately for $x$ where $A$ is a symmetric, Hermitian, skew-symmetric or skew-Hermitian linear operator and $b$ the right-hand side vector.
Usage
IterativeSolvers.minres
— Functionminres(A, b; kwargs...) -> x, [history]
Same as minres!
, but allocates a solution vector x
initialized with zeros.
IterativeSolvers.minres!
— Functionminres!(x, A, b; kwargs...) -> x, [history]
Solve Ax = b for (skew-)Hermitian matrices A using MINRES.
Arguments
x
: initial guess, will be updated in-place;A
: linear operator;b
: right-hand side.
Keywords
initially_zero::Bool = false
: iftrue
assumes thatiszero(x)
so that one matrix-vector product can be saved when computing the initial residual vector;skew_hermitian::Bool = false
: iftrue
assumes thatA
is skew-symmetric or skew-Hermitian;abstol::Real = zero(real(eltype(b)))
,reltol::Real = sqrt(eps(real(eltype(b))))
: absolute and relative tolerance for the stopping condition|r_k| ≤ max(reltol * |r_0|, abstol)
, wherer_k = A * x_k - b
is the residual in thek
th iterationNote The residual is computed only approximately.
maxiter::Int = size(A, 2)
: maximum number of iterations;log::Bool = false
: keep track of the residual norm in each iteration;verbose::Bool = false
: print convergence information during the iterations.
Return values
if log
is false
x
: approximate solution.
if log
is true
x
: approximate solution;history
: convergence history.
Implementation details
MINRES exploits the tridiagonal structure of the Hessenberg matrix. Although MINRES is mathematically equivalent to GMRES, it might not be equivalent in finite precision. MINRES updates the solution as
\[x := x_0 + (V R^{-1}) (Q^*\|r_0\|e_1)\]
where $V$ is the orthonormal basis for the Krylov subspace and $QR$ is the QR-decomposition of the Hessenberg matrix. Note that the brackets are placed slightly differently from how GMRES would update the residual.
MINRES computes $V$ and $W = VR^{-1}$ via a three-term recurrence, using only the last column of $R.$ Therefore we pre-allocate only six vectors, save only the last two entries of $Q^*\|r_0\|e_1$ and part of the last column of the Hessenberg matrix.
If $A$ is Hermitian, then the Hessenberg matrix will be real. This is exploited in the current implementation.
If $A$ is skew-Hermitian, the diagonal of the Hessenberg matrix will be imaginary, and hence we use complex arithmetic in that case.
MINRES can be used as an iterator.