API
Linear systems
IterativeRefinement.rfldiv — Functionrfldiv(A,b,f=lu; kwargs...) -> x,bnorm,bcomp,flagsCompute an accurate solution to a linear system $A x = b$ using extra-precise iterative refinement, with error bounds.
Returns solution x, a normwise relative forward error estimate bnorm, and maximum componentwise relative error estimate bcomp. Specifically, bnorm is an estimate of $‖xtrue - x‖ / ‖x‖$ (max norms). If the problem is so ill-conditioned that a good solution is unrealizable, bnorm and bcomp are set to unity (unless expert). flags contains convergence diagnostics potentially interesting to specialists.
Arguments
A: a matrix,b: a vector with the sameeltype,f: a factorization function such aslu.
Keywords
DT: higher-precision type for refinement; defaults towiden(eltype(A))verbosity: 0(default): quiet, 1: report on iterations, 2: details.equilibrate::Bool: whether the function should equilibrateA(defaulttrue).maxiter: default 20.tol: relative tolerance for convergence, in units ofeps(T).expert::Bool: whether to return questionable bounds in extreme cases.κ: the (max-norm) condition ofA(see below).F: a factorization ofA(see below).
If A has already been equilibrated, and a Factorization object F and condition estimate κ have already been computed, they may be provided as keyword arguments; no check for consistency is done here.
Uses the algorithm of Demmel et al. ACM TOMS, 32, 325 (2006).
IterativeRefinement.equilibrators — Functionequilibrators(A) -> R,Ccompute row- and column-wise scaling vectors R,C for a matrix A such that the absolute value of the largest element in any row or column of Diagonal(R)*A*Diagonal(C) is close to unity. Designed to reduce the condition number of the working matrix.
IterativeRefinement.condInfest — FunctioncondInfest(A,F,anorm)computes an approximation to the condition of matrix A in the infinity-norm, using factorization F and the precomputed infinity norm anorm of A.
Eigensystems
IterativeRefinement.rfeigen — Functionrfeigen(A,x,λ,DT) => λnew, xnew, statusImprove the precision of a computed eigenpair (x,λ) for matrix A via multi-precision iterative refinement, using more-precise real type DT.
The higher precision DT is only used for residual computation (i.e. matrix-vector products), so this can be much faster than a full eigensystem solution with precise eltype. This method works on a single eigenpair, and can fail spectacularly if there is another eigenvalue nearby.
rfeigen(A,λ,DT) => λnew, xnew, statusLike rfeigen(A,x,λ,DT), but initialize x via one step of inverse iteration.
rfeigen(A, S::Schur, idxλ, DT, maxiter=5) -> vals, vecs, statusImproves the precision of a cluster of eigenvalues of square matrix A via multi-precision iterative refinement, using more-precise real type DT, using a pre-computed Schur decomposition S. Returns improved estimates of eigenvalues and vectors generating the corresponding invariant subspace.
This method works on the set of eigenvalues in S.values indexed by idxλ. It is designed to handle (nearly) defective cases, but will fail if the matrix is extremely non-normal or the initial estimates are poor.
If S is a quasi-triangular "real Schur", it will be converted to a complex upper-triangular Schur decomposition. S may be a partial decomposition (i.e. fewer than size(A,1) vectors).