Arnoldi Iteration

ExponentialUtilities.arnoldiFunction
arnoldi(A,b[;m,tol,opnorm,iop]) -> Ks

Performs m arnoldi iterations to obtain the Krylov subspace K_m(A,b).

The n x (m + 1) basis vectors getV(Ks) and the (m + 1) x m upper Hessenberg matrix getH(Ks) are related by the recurrence formula

v_1=b,\quad Av_j = \sum_{i=1}^{j+1}h_{ij}v_i\quad(j = 1,2,\ldots,m)

iop determines the length of the incomplete orthogonalization procedure [1]. The default value of 0 indicates full Arnoldi. For symmetric/Hermitian A, iop will be ignored and the Lanczos algorithm will be used instead.

Refer to KrylovSubspace for more information regarding the output.

Happy-breakdown occurs whenever norm(v_j) < tol * opnorm, in this case, the dimension of Ks is smaller than m.

source

API

ExponentialUtilities.KrylovSubspaceType
KrylovSubspace{T}(n,[maxiter=30]) -> Ks
KrylovSubspace{T, U, VType}(n,[maxiter=30]) -> Ks

Constructs an uninitialized Krylov subspace, which can be filled by arnoldi!.

The dimension of the subspace, Ks.m, can be dynamically altered but should be smaller than maxiter, the maximum allowed arnoldi iterations.

The type of the (extended) orthonormal basis vector matrix V may be specified as VType. This is required e. g. for GPUArrays.

U determines eltype(H).

getV(Ks) -> V
getH(Ks) -> H

Access methods for the (extended) orthonormal basis V and the (extended) Gram-Schmidt coefficients H. Both methods return a view into the storage arrays and has the correct dimensions as indicated by Ks.m.

resize!(Ks, maxiter) -> Ks

Resize Ks to a different maxiter, destroying its contents.

This is an expensive operation and should be used scarcely.

source
ExponentialUtilities.firststep!Function
firststep!(Ks, V, H, b) -> nothing

Compute the first step of Arnoldi or Lanczos iteration.

source
firststep!(Ks, V, H, b, b_aug, t, mu, l) -> nothing

Compute the first step of Arnoldi or Lanczos iteration of augmented system.

source
ExponentialUtilities.coeffFunction
coeff(::Type,α)

Helper functions that returns the real part if that is what is required (for Hermitian matrices), otherwise returns the value as-is.

source
  • 1Koskela, A. (2015). Approximating the matrix exponential of an advection-diffusion operator using the incomplete orthogonalization method. In Numerical Mathematics and Advanced Applications-ENUMATH 2013 (pp. 345-353). Springer, Cham.