arnoldi Interface

public interface arnoldi

Description

Given a square linear operator , find matrices and such that

where is an orthogonal basis and is upper Hessenberg.

Algorithmic Features

  • The operator only needs to be accessed through matrix-vector products.
  • Constructs an orthonormal Krylov basis via the Gram-Schmidt process.
  • Constructs an upper Hessenberg matrix whose eigenvalues approximates those of .
  • Checks for convergence and invariant subspaces.

References

  • Y. Saad. "Iterative methods for sparse linear systems", SIAM 2nd edition, 2003. see Chapter 6.3: Arnoldi's method.

Syntax

    call arnoldi(A, X, H, info [, kstart] [, kend] [, tol] [, transpose] [, blksize])

Arguments

A : Linear operator derived from one the base types provided by the AbstractLinops module. The operator needs to be square, i.e. the dimension of its domain and co-domain is the same. It is an intent(in) argument.

X : Array of types derived from one the base types provided by the AbstractVectors module. It needs to be consistent with the type of A. On exit, it contains the the computed Krylov vectors. The first entry X(1) is the starting vector for the Arnoldi factorization. Additionally, the maximum number of Arnoldi steps is equal to size(X) - 1. It is an intent(inout) argument.

H : real or complex rank-2 array. On exit, it contains the upper Hessenberg matrix computed from the Arnoldi factorization. It is an intent(inout) argument.

info : integer variable. It is the LightKrylov information flag. On exit, if info > 0, the Arnoldi factorization experienced a lucky breakdown. The array of Krylov vectors X spans an -invariant subpsace of dimension info.

kstart (optional): integer value determining the index of the first Arnoldi step to be computed. By default, kstart = 1.

kend (optional): integer value determining the index of the last Arnoldi step to be computed. By default, kend = size(X) - 1.

tol (optional): Numerical tolerance below which a subspace is considered to be -invariant. By default tol = atol_sp or tol = atol_rp depending on the kind of A.

transpose (optional): logical flag determining whether the Arnoldi factorization is applied to or . Default transpose = .false.

blksize (optional): integer value determining the dimension of a block for the block Arnoldi factorization. Default is blksize=1.


Module Procedures

private subroutine arnoldi_rsp(A, X, H, info, kstart, kend, tol, transpose, blksize)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_rsp), intent(in) :: A

Linear operator to be factorized.

class(abstract_vector_rsp), intent(inout) :: X(:)

Orthogonal basis for the generated Krylov subspace.

real(kind=sp), intent(inout) :: H(:,:)

Upper Hessenberg matrix.

integer, intent(out) :: info

Information flag.

integer, intent(in), optional :: kstart

Starting index for the Arnoldi factorization (default 1).

integer, intent(in), optional :: kend

Final index for the Arnoldi factorization (default size(X)-1)

real(kind=sp), intent(in), optional :: tol

Tolerance to determine whether an invariant subspace has been computed or not.

logical, intent(in), optional :: transpose

Whether is being transposed or not (default .false.)

integer, intent(in), optional :: blksize

Block size for block Arnoldi (default 1).

private subroutine arnoldi_rdp(A, X, H, info, kstart, kend, tol, transpose, blksize)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_rdp), intent(in) :: A

Linear operator to be factorized.

class(abstract_vector_rdp), intent(inout) :: X(:)

Orthogonal basis for the generated Krylov subspace.

real(kind=dp), intent(inout) :: H(:,:)

Upper Hessenberg matrix.

integer, intent(out) :: info

Information flag.

integer, intent(in), optional :: kstart

Starting index for the Arnoldi factorization (default 1).

integer, intent(in), optional :: kend

Final index for the Arnoldi factorization (default size(X)-1)

real(kind=dp), intent(in), optional :: tol

Tolerance to determine whether an invariant subspace has been computed or not.

logical, intent(in), optional :: transpose

Whether is being transposed or not (default .false.)

integer, intent(in), optional :: blksize

Block size for block Arnoldi (default 1).

private subroutine arnoldi_csp(A, X, H, info, kstart, kend, tol, transpose, blksize)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_csp), intent(in) :: A

Linear operator to be factorized.

class(abstract_vector_csp), intent(inout) :: X(:)

Orthogonal basis for the generated Krylov subspace.

complex(kind=sp), intent(inout) :: H(:,:)

Upper Hessenberg matrix.

integer, intent(out) :: info

Information flag.

integer, intent(in), optional :: kstart

Starting index for the Arnoldi factorization (default 1).

integer, intent(in), optional :: kend

Final index for the Arnoldi factorization (default size(X)-1)

real(kind=sp), intent(in), optional :: tol

Tolerance to determine whether an invariant subspace has been computed or not.

logical, intent(in), optional :: transpose

Whether is being transposed or not (default .false.)

integer, intent(in), optional :: blksize

Block size for block Arnoldi (default 1).

private subroutine arnoldi_cdp(A, X, H, info, kstart, kend, tol, transpose, blksize)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_cdp), intent(in) :: A

Linear operator to be factorized.

class(abstract_vector_cdp), intent(inout) :: X(:)

Orthogonal basis for the generated Krylov subspace.

complex(kind=dp), intent(inout) :: H(:,:)

Upper Hessenberg matrix.

integer, intent(out) :: info

Information flag.

integer, intent(in), optional :: kstart

Starting index for the Arnoldi factorization (default 1).

integer, intent(in), optional :: kend

Final index for the Arnoldi factorization (default size(X)-1)

real(kind=dp), intent(in), optional :: tol

Tolerance to determine whether an invariant subspace has been computed or not.

logical, intent(in), optional :: transpose

Whether is being transposed or not (default .false.)

integer, intent(in), optional :: blksize

Block size for block Arnoldi (default 1).