LightKrylov_IterativeSolvers Module

This module provides some of the most important computational routines provided by LightKrylov. These include:

  • eigs : Compute the leading eigenpairs of a square linear operator .
  • eighs : Compute the leading eigenpairs of a symmetric positive definite operator .
  • svds : Compute the leading singular triplets of a linear operator .
  • gmres : Solve the linear system using the generalized minimum residual method.
  • cg : Solve the linear system where is symmetric positive definite using the Conjugate Gradient method.

It also provides abstract interfaces to pass user-defined solvers and preconditioners to LightKrylov. Note that these features are still experimental however.



Interfaces

public interface cg

Description

Given a symmetric (positive definite) matrix , solves the linear system

using the Conjugate Gradient method.

References

  • Hestenes, M. R., and Stiefel, E. (1952). "Methods of Conjugate Gradients for Solving Linear Systems," Journal of Research of the National Bureau of Standards, 49(6), 409–436.

Syntax

    call cg(A, b, x, info [, rtol] [, atol] [, preconditioner] [, options])

Arguments

  • A : Linear operator derived from one of the abstract_sym_linop or abstract_hermitian_linop types provided

  • b : Right-hand side vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(in) argument.

  • x : On entry, initial guess for the solution. On exit, the solution computed by cg. It is a vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(inout) argument.

  • info : integer information flag.

  • rtol (optional) : real relative tolerance for the solver.

  • atol (optional) : real absolute tolerance for the solver.

  • preconditioner (optional) : Right preconditioner used to solve the system. It needs to be consistent with the abstract_preconditioner interface. It is an intent(in) argument.

  • options (optional) : Container for the gmres options given by the cg_opts type. It is an intent(in) argument.

  • private module subroutine cg_cdp(A, b, x, info, rtol, atol, preconditioner, options, meta)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_cdp), intent(inout) :: A

    Linear operator to be inverted.

    class(abstract_vector_cdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner (not yet supported).

    type(cg_dp_opts), intent(in), optional :: options

    Options for the conjugate gradient solver.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine cg_csp(A, b, x, info, rtol, atol, preconditioner, options, meta)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_csp), intent(inout) :: A

    Linear operator to be inverted.

    class(abstract_vector_csp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner (not yet supported).

    type(cg_sp_opts), intent(in), optional :: options

    Options for the conjugate gradient solver.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine cg_rdp(A, b, x, info, rtol, atol, preconditioner, options, meta)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rdp), intent(inout) :: A

    Linear operator to be inverted.

    class(abstract_vector_rdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner (not yet supported).

    type(cg_dp_opts), intent(in), optional :: options

    Options for the conjugate gradient solver.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine cg_rsp(A, b, x, info, rtol, atol, preconditioner, options, meta)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rsp), intent(inout) :: A

    Linear operator to be inverted.

    class(abstract_vector_rsp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner (not yet supported).

    type(cg_sp_opts), intent(in), optional :: options

    Options for the conjugate gradient solver.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

public interface eighs

Description

Computes the leading eigenpairs of a symmetric operator using the Lanczos iterative process. Given a square linear operator , it finds the leading eigvalues and eigvectors such that:

The subspace is constructed via Lanczos factorization, resulting in a symmetric tridiagonal matrix . The eigenvalues of are approximated by those of and the eigenvectors are computed accordingly.

References

  • Lanczos, C. (1950). "An Iteration Method for the Solution of the Eigenvalue Problem of Linear Differential and Integral Operators". United States Governm. Press Office.

Syntax

    call eighs(A, X, eigvals, residuals, info [, kdim] [,tolerance])

Arguments

  • A : Linear operator derived from abstract_sym_linop_rsp, abstract_sym_linop_rdp, abstract_hermitian_linop_csp or abstract_hermitian_linop_cdp whose leading eigenpairs need to be computed. It is an intent(inout) argument.

  • X : Array of abstract_vectors with the same type and kind as A. On exit, it contains the leading eigenvectors of A. Note that the dimension of X fixes the number of eigenpairs computed.

  • eigvals : Rank-1 array of real numbers. On exit, it contains the leading eigenvalues of A. It is an intent(out) argument.

  • residuals : Rank-1 array of real numbers. On exit, it contains the residuals associated with each eigenpairs. It is an intent(out) argument.

  • info : integer Information flag.

  • kdim (optional) : integer, maximum dimension of the Krylov subspace used to approximate the leading eigenpairs. It is an intent(in) argument. By default, kdim = 4*size(X).

  • tolerance (optional) : real tolerance below which an eigenpair is considered as being converged. It is an intent(in) agument. By default, tolerance = rtol_sp or tolerance = rtol_dp.

Note

This implementation does not currently include an automatic restarting procedure such as krylov_schur for eigs. This is work in progress.

  • private module subroutine eighs_cdp(A, X, eigvals, residuals, info, x0, kdim, tolerance, write_intermediate)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_cdp), intent(inout) :: A

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigevectors of .

    real(kind=dp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpairs.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_cdp), intent(in), optional :: x0

    Optional starting vector to generate the Krylov subspace.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine eighs_csp(A, X, eigvals, residuals, info, x0, kdim, tolerance, write_intermediate)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_csp), intent(inout) :: A

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigevectors of .

    real(kind=sp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpairs.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_csp), intent(in), optional :: x0

    Optional starting vector to generate the Krylov subspace.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine eighs_rdp(A, X, eigvals, residuals, info, x0, kdim, tolerance, write_intermediate)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rdp), intent(inout) :: A

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigevectors of .

    real(kind=dp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpairs.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rdp), intent(in), optional :: x0

    Optional starting vector to generate the Krylov subspace.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine eighs_rsp(A, X, eigvals, residuals, info, x0, kdim, tolerance, write_intermediate)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rsp), intent(inout) :: A

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigevectors of .

    real(kind=sp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpairs.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rsp), intent(in), optional :: x0

    Optional starting vector to generate the Krylov subspace.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

public interface eigs

Description

Computes the leading eigenpairs of a square linear operator using the Arnoldi iterative process. Given a square linear operator , it finds the leading eigvalues and eigvectorss such that:

or

The subspace is constructed via Arnoldi factorization, resulting in an upper Hessenberg matrix . The eigenvalues of are approximated by those of and the eigenvectors are computed accordingly.

References

  • Arnoldi, W. E. (1951). "The Principle of Minimized Iterations in the Solution of the Matrix Eigenvalue Problem." Quarterly of Applied Mathematics, 9(1), 17–29.

Syntax

    call eigs(A, X, eigvals, residuals, info [, kdim] [, select] [,tolerance] [, transpose])

Arguments

  • A : Linear operator derived from abstract_sym_linop_rsp, abstract_sym_linop_rdp, abstract_hermitian_linop_csp or abstract_hermitian_linop_cdp whose leading eigenpairs need to be computed. It is an intent(inout) argument.

  • X : Array of abstract_vectors with the same type and kind as A. On exit, it contains the leading eigenvectors of A. Note that the dimension of X fixes the number of eigenpairs computed.

  • eigvals : Rank-1 array of real numbers. On exit, it contains the leading eigenvalues of A. It is an intent(out) argument.

  • residuals : Rank-1 array of real numbers. On exit, it contains the residuals associated with each eigenpairs. It is an intent(out) argument.

  • info : integer Information flag.

  • kdim (optional) : integer, maximum dimension of the Krylov subspace used to approximate the leading eigenpairs. It is an intent(in) argument. By default, kdim = 4*size(X).

  • select (optional) : Function to select which eigenvalues to compute.

  • tolerance (optional) : real tolerance below which an eigenpair is considered as being converged. It is an intent(in) agument. By default, tolerance = rtol_sp or tolerance = rtol_dp.

  • transpose (optional) : logical flag determining whether the eigenvalues of or need to be computed.

  • private subroutine eigs_rsp(A, X, eigvals, residuals, info, x0, kdim, tolerance, transpose, write_intermediate)

    Arguments

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

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigenvectors of .

    complex(kind=sp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rsp), intent(in), optional :: x0

    Optional starting vector for generating the Krylov subspace.

    integer, intent(in), optional :: kdim

    Maximum dimension of the Krylov subspace (optional).

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private subroutine eigs_rdp(A, X, eigvals, residuals, info, x0, kdim, tolerance, transpose, write_intermediate)

    Arguments

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

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigenvectors of .

    complex(kind=dp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rdp), intent(in), optional :: x0

    Optional starting vector for generating the Krylov subspace.

    integer, intent(in), optional :: kdim

    Maximum dimension of the Krylov subspace (optional).

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private subroutine eigs_csp(A, X, eigvals, residuals, info, x0, kdim, tolerance, transpose, write_intermediate)

    Arguments

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

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigenvectors of .

    complex(kind=sp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_csp), intent(in), optional :: x0

    Optional starting vector for generating the Krylov subspace.

    integer, intent(in), optional :: kdim

    Maximum dimension of the Krylov subspace (optional).

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private subroutine eigs_cdp(A, X, eigvals, residuals, info, x0, kdim, tolerance, transpose, write_intermediate)

    Arguments

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

    Linear operator whose leading eigenpairs need to be computed.

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

    Leading eigenvectors of .

    complex(kind=dp), intent(out), allocatable :: eigvals(:)

    Leading eigenvalues of .

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_cdp), intent(in), optional :: x0

    Optional starting vector for generating the Krylov subspace.

    integer, intent(in), optional :: kdim

    Maximum dimension of the Krylov subspace (optional).

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

public interface fgmres

Description

Solve a square linear system of equations

using the Flexible Generalized Minimum RESidual (FGMRES) method.

References

  • Saad Y. and Schultz M. H. "GMRES: A generalized minimal residual algorithm for solving nonsymmetric linear systems." SIAM Journal on Scientific and Statistical Computing, 7(3), 1986.

Syntax

    call fgmres(A, b, x, info [, rtol] [, atol] [, preconditioner] [, options] [, transpose])

Arguments

  • A : Linear operator derived from one of the abstract_linop provided by the AbstractLinops module. It is an intent(inout) argument.

  • b : Right-hand side vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(in) argument.

  • x : On entry, initial guess for the solution. On exit, the solution computed by gmres. It is a vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(inout) argument.

  • info : integer information flag.

  • rtol (optional) : real relative tolerance for the solver.

  • atol (optional) : real absolute tolerance for the solver.

  • preconditioner (optional) : Right preconditioner used to solve the system. It needs to be consistent with the abstract_preconditioner interface. It is an intent(in) argument.

  • options (optional) : Container for the gmres options given by the gmres_opts type. It is an intent(in) argument.

  • transpose (optional): logical flag controlling whether or is being solved.

  • private module subroutine dense_fgmres_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=dp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    complex(kind=dp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_fgmres_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=sp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    complex(kind=sp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_fgmres_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    real(kind=dp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_fgmres_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=sp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    real(kind=sp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine fgmres_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_cdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine fgmres_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_csp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine fgmres_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_rdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine fgmres_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_rsp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

public interface gmres

Description

Solve a square linear system of equations

using the Generalized Minimum RESidual (GMRES) method.

References

  • Saad Y. and Schultz M. H. "GMRES: A generalized minimal residual algorithm for solving nonsymmetric linear systems." SIAM Journal on Scientific and Statistical Computing, 7(3), 1986.

Syntax

    call gmres(A, b, x, info [, rtol] [, atol] [, preconditioner] [, options] [, transpose])

Arguments

  • A : Linear operator derived from one of the abstract_linop provided by the AbstractLinops module. It is an intent(inout) argument.

  • b : Right-hand side vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(in) argument.

  • x : On entry, initial guess for the solution. On exit, the solution computed by gmres. It is a vector derived from one the abstract_vector types provided by the AbstractVectors module. It needs to have the same type and kind as A. It is an intent(inout) argument.

  • info : integer information flag.

  • rtol (optional) : real relative tolerance for the solver.

  • atol (optional) : real absolute tolerance for the solver.

  • preconditioner (optional) : Right preconditioner used to solve the system. It needs to be consistent with the abstract_preconditioner interface. It is an intent(in) argument.

  • options (optional) : Container for the gmres options given by the gmres_opts type. It is an intent(in) argument.

  • transpose (optional): logical flag controlling whether or is being solver.

  • meta (optional) : Container for the gmres metada. It needs to be of type gmres_medata.

  • private module subroutine dense_gmres_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=dp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    complex(kind=dp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_gmres_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=sp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    complex(kind=sp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_gmres_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    real(kind=dp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine dense_gmres_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=sp), intent(in) :: A(:,:)

    Linear operator to be inverted.

    real(kind=sp), intent(in) :: b(:)

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine gmres_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_cdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine gmres_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_csp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine gmres_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_rdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

  • private module subroutine gmres_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Arguments

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

    Linear operator to be inverted.

    class(abstract_vector_rsp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner (optional).

    class(abstract_opts), intent(in), optional :: options

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

public interface save_eigenspectrum

Description

Utility function to save the eigenspectrum computed from the Arnoldi factorization. It outpost a .npy file.

Syntax

    call save_eigenspectrum(eigvals, residuals, fname)

Arguments

  • eigvals : complex rank-1 array containing the eigenvalues.

  • residuals : real rank-1 array containing the residuals associated to each eigenvalues.

fname : Name of the file to save the eigenspectrum.

  • private module subroutine save_eigenspectrum_cdp(lambda, residuals, fname)

    Saves the eigenspectrum and corresponding residuals to disk use the npy binary format.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=dp), intent(in) :: lambda(:)

    Eigenalues.

    real(kind=dp), intent(in) :: residuals(:)

    Residual of the corresponding Ritz eigenpairs.

    character(len=*), intent(in) :: fname

    Name of the output file.

  • private module subroutine save_eigenspectrum_csp(lambda, residuals, fname)

    Saves the eigenspectrum and corresponding residuals to disk use the npy binary format.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=sp), intent(in) :: lambda(:)

    Eigenalues.

    real(kind=sp), intent(in) :: residuals(:)

    Residual of the corresponding Ritz eigenpairs.

    character(len=*), intent(in) :: fname

    Name of the output file.

  • private module subroutine save_eigenspectrum_rdp(lambda, residuals, fname)

    Saves the eigenspectrum and corresponding residuals to disk use the npy binary format.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: lambda(:)

    Eigenalues.

    real(kind=dp), intent(in) :: residuals(:)

    Residual of the corresponding Ritz eigenpairs.

    character(len=*), intent(in) :: fname

    Name of the output file.

  • private module subroutine save_eigenspectrum_rsp(lambda, residuals, fname)

    Saves the eigenspectrum and corresponding residuals to disk use the npy binary format.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=sp), intent(in) :: lambda(:)

    Eigenalues.

    real(kind=sp), intent(in) :: residuals(:)

    Residual of the corresponding Ritz eigenpairs.

    character(len=*), intent(in) :: fname

    Name of the output file.

public interface svds

Description

Computes the leading singular triplets of an arbitrary linear operator using the Lanczos iterative process. Given a linear operator , it finds the leading singular values and singular vectors such that:

The subspaces and are constructed via Lanczos factorization, resulting in a bidiagonal matrix . The singular values of are approximated by those of and the singular vectors are computed accordingly.

References

  • Golub, G. H., & Kahan, W. (1965). "Calculating the Singular Values and Pseudo-Inverse of a Matrix."
  • Baglama, J., & Reichel, L. (2005). "Augmented implicitly restarted Lanczos bidiagonalization methods."
  • R. M. Larsen. "Lanczos bidiagonalization with partial reorthogonalization." Technical Report, 1998.

Syntax

    call svds(A, U, S, V, residuals, info [, kdim] [,tolerance])

Arguments

  • A : Linear operator derived from abstract_sym_linop_rsp, abstract_sym_linop_rdp, abstract_hermitian_linop_csp or abstract_hermitian_linop_cdp whose leading eigenpairs need to be computed. It is an intent(inout) argument.

  • U : Array of abstract_vectors with the same type and kind as A. On exit, it contains the left singular vectors of A. Note that the dimension of U fixes the number of eigenpairs computed.

  • S : Rank-1 array of real numbers. On exit, it contains the leading singular values of A. It is an intent(out) argument.

  • V : Array of abstract_vectors with the same type and kind as A. On exit, it contains the left singular vectors of A. Note that the dimension of U fixes the number of eigenpairs computed.

  • residuals : Rank-1 array of real numbers. On exit, it contains the residuals associated with each singular triplet. It is an intent(out) argument.

  • info : integer Information flag.

  • kdim (optional) : integer, maximum dimension of the Krylov subspace used to approximate the leading singular triplets. It is an intent(in) argument. By default, kdim = 4*size(X).

  • tolerance (optional) : real tolerance below which a triplet is considered as being converged. It is an intent(in) agument. By default, tolerance = rtol_sportolerance = rtol_dp`.

Note

This implementation does not currently include an automatic restarting procedure such as krylov_schur for eigs. This is work in progress.

  • private module subroutine svds_cdp(A, U, S, V, residuals, info, u0, kdim, tolerance, write_intermediate)

    Arguments

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

    Linear operator whose leading singular triplets need to be computed.

    class(abstract_vector_cdp), intent(out) :: U(:)

    Leading left singular vectors.

    real(kind=dp), intent(out), allocatable :: S(:)

    Leading singular values.

    class(abstract_vector_cdp), intent(out) :: V(:)

    Leading right singular vectors.

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_cdp), intent(in), optional :: u0
    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine svds_csp(A, U, S, V, residuals, info, u0, kdim, tolerance, write_intermediate)

    Arguments

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

    Linear operator whose leading singular triplets need to be computed.

    class(abstract_vector_csp), intent(out) :: U(:)

    Leading left singular vectors.

    real(kind=sp), intent(out), allocatable :: S(:)

    Leading singular values.

    class(abstract_vector_csp), intent(out) :: V(:)

    Leading right singular vectors.

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_csp), intent(in), optional :: u0
    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine svds_rdp(A, U, S, V, residuals, info, u0, kdim, tolerance, write_intermediate)

    Arguments

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

    Linear operator whose leading singular triplets need to be computed.

    class(abstract_vector_rdp), intent(out) :: U(:)

    Leading left singular vectors.

    real(kind=dp), intent(out), allocatable :: S(:)

    Leading singular values.

    class(abstract_vector_rdp), intent(out) :: V(:)

    Leading right singular vectors.

    real(kind=dp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rdp), intent(in), optional :: u0
    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?

  • private module subroutine svds_rsp(A, U, S, V, residuals, info, u0, kdim, tolerance, write_intermediate)

    Arguments

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

    Linear operator whose leading singular triplets need to be computed.

    class(abstract_vector_rsp), intent(out) :: U(:)

    Leading left singular vectors.

    real(kind=sp), intent(out), allocatable :: S(:)

    Leading singular values.

    class(abstract_vector_rsp), intent(out) :: V(:)

    Leading right singular vectors.

    real(kind=sp), intent(out), allocatable :: residuals(:)

    Residuals associated to each Ritz eigenpair.

    integer, intent(out) :: info

    Information flag.

    class(abstract_vector_rsp), intent(in), optional :: u0
    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: write_intermediate

    Write intermediate eigenvalues to file during iteration?


Abstract Interfaces

abstract interface

  • public subroutine abstract_linear_solver_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Abstract interface to use a user-defined linear solver in LightKrylov.

    Arguments

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

    Linear operator to invert.

    class(abstract_vector_cdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag. In case of successful exit, the flag should return the number of iterations required for convergence.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_cdp), intent(inout), optional :: preconditioner

    Preconditioner.

    class(abstract_opts), intent(in), optional :: options

    Options passed to the linear solver.

    logical, intent(in), optional :: transpose

    Determine whether (.false.) or (.true.) is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

abstract interface

  • public subroutine abstract_linear_solver_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Abstract interface to use a user-defined linear solver in LightKrylov.

    Arguments

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

    Linear operator to invert.

    class(abstract_vector_csp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag. In case of successful exit, the flag should return the number of iterations required for convergence.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_csp), intent(inout), optional :: preconditioner

    Preconditioner.

    class(abstract_opts), intent(in), optional :: options

    Options passed to the linear solver.

    logical, intent(in), optional :: transpose

    Determine whether (.false.) or (.true.) is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

abstract interface

  • public subroutine abstract_linear_solver_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Abstract interface to use a user-defined linear solver in LightKrylov.

    Arguments

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

    Linear operator to invert.

    class(abstract_vector_rdp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag. In case of successful exit, the flag should return the number of iterations required for convergence.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rdp), intent(inout), optional :: preconditioner

    Preconditioner.

    class(abstract_opts), intent(in), optional :: options

    Options passed to the linear solver.

    logical, intent(in), optional :: transpose

    Determine whether (.false.) or (.true.) is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.

abstract interface

  • public subroutine abstract_linear_solver_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose, meta)

    Abstract interface to use a user-defined linear solver in LightKrylov.

    Arguments

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

    Linear operator to invert.

    class(abstract_vector_rsp), intent(in) :: b

    Right-hand side vector.

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

    Solution vector.

    integer, intent(out) :: info

    Information flag. In case of successful exit, the flag should return the number of iterations required for convergence.

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

    Relative solver tolerance

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

    Absolute solver tolerance

    class(abstract_precond_rsp), intent(inout), optional :: preconditioner

    Preconditioner.

    class(abstract_opts), intent(in), optional :: options

    Options passed to the linear solver.

    logical, intent(in), optional :: transpose

    Determine whether (.false.) or (.true.) is being used.

    class(abstract_metadata), intent(out), optional :: meta

    Metadata.


Derived Types

type, public, abstract ::  abstract_precond_cdp

Type-Bound Procedures

procedure(abstract_apply_cdp), public, deferred, pass(self) :: apply

type, public, abstract ::  abstract_precond_csp

Type-Bound Procedures

procedure(abstract_apply_csp), public, deferred, pass(self) :: apply

type, public, abstract ::  abstract_precond_rdp

Type-Bound Procedures

procedure(abstract_apply_rdp), public, deferred, pass(self) :: apply

type, public, abstract ::  abstract_precond_rsp

Type-Bound Procedures

procedure(abstract_apply_rsp), public, deferred, pass(self) :: apply

type, public, extends(abstract_metadata) ::  cg_dp_metadata

Conjugate gradient metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag

integer, public :: info = 0

Copy of the information flag for completeness

integer, public :: n_iter = 0

Iteration counter

real(kind=dp), public, dimension(:), allocatable :: res

Residual history

Type-Bound Procedures

procedure, public, pass(self) :: print => print_cg_dp
procedure, public, pass(self) :: reset => reset_cg_dp

type, public, extends(abstract_opts) ::  cg_dp_opts

Conjugate gradient options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print interation metadata on exit (default = .false.)

integer, public :: maxiter = 100

Maximum number of cg iterations (default: 100).

type, public, extends(abstract_metadata) ::  cg_sp_metadata

Conjugate gradient metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag

integer, public :: info = 0

Copy of the information flag for completeness

integer, public :: n_iter = 0

Iteration counter

real(kind=sp), public, dimension(:), allocatable :: res

Residual history

Type-Bound Procedures

procedure, public, pass(self) :: print => print_cg_sp
procedure, public, pass(self) :: reset => reset_cg_sp

type, public, extends(abstract_opts) ::  cg_sp_opts

Conjugate gradient options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print interation metadata on exit (default = .false.)

integer, public :: maxiter = 100

Maximum number of cg iterations (default: 100).

type, public, extends(abstract_metadata) ::  fgmres_dp_metadata

FGMRES metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag.

integer, public :: info = 0

Copy of the information flag for completeness.

integer, public :: n_inner = 0

Number of inner iterations.

integer, public :: n_iter = 0

Total iteration counter.

integer, public :: n_outer = 0

Number of outer iterations (i.e. restarts)

real(kind=dp), public, dimension(:), allocatable :: res

Residual history.

Type-Bound Procedures

procedure, public, pass(self) :: print => print_fgmres_dp
procedure, public, pass(self) :: reset => reset_fgmres_dp

type, public, extends(abstract_opts) ::  fgmres_dp_opts

FGMRES options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print iteration metadata on exit (default: .false.).

integer, public :: kdim = 30

Dimension of the Krylov subspace (default: 30).

integer, public :: maxiter = 10

Maximum number of fgmres restarts (default: 10).

logical, public :: sanity_check = .true.

Performs extra matrix-vector product for sanity check.

type, public, extends(abstract_metadata) ::  fgmres_sp_metadata

FGMRES metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag.

integer, public :: info = 0

Copy of the information flag for completeness.

integer, public :: n_inner = 0

Number of inner iterations.

integer, public :: n_iter = 0

Total iteration counter.

integer, public :: n_outer = 0

Number of outer iterations (i.e. restarts)

real(kind=sp), public, dimension(:), allocatable :: res

Residual history.

Type-Bound Procedures

procedure, public, pass(self) :: print => print_fgmres_sp
procedure, public, pass(self) :: reset => reset_fgmres_sp

type, public, extends(abstract_opts) ::  fgmres_sp_opts

FGMRES options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print iteration metadata on exit (default: .false.).

integer, public :: kdim = 30

Dimension of the Krylov subspace (default: 30).

integer, public :: maxiter = 10

Maximum number of fgmres restarts (default: 10).

logical, public :: sanity_check = .true.

Performs extra matrix-vector product for sanity check.

type, public, extends(abstract_metadata) ::  gmres_dp_metadata

GMRES metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag.

integer, public :: info = 0

Copy of the information flag for completeness.

integer, public :: n_inner = 0

Number of inner iterations.

integer, public :: n_iter = 0

Total iteration counter.

integer, public :: n_outer = 0

Number of outer iterations (i.e. restarts)

real(kind=dp), public, dimension(:), allocatable :: res

Residual history.

Type-Bound Procedures

procedure, public, pass(self) :: print => print_gmres_dp
procedure, public, pass(self) :: reset => reset_gmres_dp

type, public, extends(abstract_opts) ::  gmres_dp_opts

GMRES options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print iteration metadata on exit (default: .false.).

integer, public :: kdim = 30

Dimension of the Krylov subspace (default: 30).

integer, public :: maxiter = 10

Maximum number of gmres restarts (default: 10).

logical, public :: sanity_check = .true.

Performs extra matrix-vector product for sanity check.

type, public, extends(abstract_metadata) ::  gmres_sp_metadata

GMRES metadata.

Components

Type Visibility Attributes Name Initial
logical, public :: converged = .false.

Convergence flag.

integer, public :: info = 0

Copy of the information flag for completeness.

integer, public :: n_inner = 0

Number of inner iterations.

integer, public :: n_iter = 0

Total iteration counter.

integer, public :: n_outer = 0

Number of outer iterations (i.e. restarts)

real(kind=sp), public, dimension(:), allocatable :: res

Residual history.

Type-Bound Procedures

procedure, public, pass(self) :: print => print_gmres_sp
procedure, public, pass(self) :: reset => reset_gmres_sp

type, public, extends(abstract_opts) ::  gmres_sp_opts

GMRES options.

Components

Type Visibility Attributes Name Initial
logical, public :: if_print_metadata = .false.

Print iteration metadata on exit (default: .false.).

integer, public :: kdim = 30

Dimension of the Krylov subspace (default: 30).

integer, public :: maxiter = 10

Maximum number of gmres restarts (default: 10).

logical, public :: sanity_check = .true.

Performs extra matrix-vector product for sanity check.


Subroutines

public subroutine write_results_cdp(filename, vals, res, tol)

Prints the intermediate results of iterative eigenvalue/singular value decompositions

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename

Output filename. This file will be overwritten

complex(kind=dp), intent(in) :: vals(:)

Intermediate values

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

Residuals

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

Convergence tolerance

public subroutine write_results_csp(filename, vals, res, tol)

Prints the intermediate results of iterative eigenvalue/singular value decompositions

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename

Output filename. This file will be overwritten

complex(kind=sp), intent(in) :: vals(:)

Intermediate values

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

Residuals

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

Convergence tolerance

public subroutine write_results_rdp(filename, vals, res, tol)

Prints the intermediate results of iterative eigenvalue/singular value decompositions

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename

Output filename. This file will be overwritten

real(kind=dp), intent(in) :: vals(:)

Intermediate values

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

Residuals

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

Convergence tolerance

public subroutine write_results_rsp(filename, vals, res, tol)

Prints the intermediate results of iterative eigenvalue/singular value decompositions

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename

Output filename. This file will be overwritten

real(kind=sp), intent(in) :: vals(:)

Intermediate values

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

Residuals

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

Convergence tolerance


Module Procedures

module procedure /github/workspace/API-doc/module/lightkrylov_iterativesolvers.html save_eigenspectrum_cdp private module subroutine save_eigenspectrum_cdp(lambda, residuals, fname)

Arguments

Type IntentOptional Attributes Name
complex(kind=dp), intent(in) :: lambda(:)

Eigenalues.

real(kind=dp), intent(in) :: residuals(:)

Residual of the corresponding Ritz eigenpairs.

character(len=*), intent(in) :: fname

Name of the output file.

module procedure /github/workspace/API-doc/module/lightkrylov_iterativesolvers.html save_eigenspectrum_csp private module subroutine save_eigenspectrum_csp(lambda, residuals, fname)

Arguments

Type IntentOptional Attributes Name
complex(kind=sp), intent(in) :: lambda(:)

Eigenalues.

real(kind=sp), intent(in) :: residuals(:)

Residual of the corresponding Ritz eigenpairs.

character(len=*), intent(in) :: fname

Name of the output file.

module procedure /github/workspace/API-doc/module/lightkrylov_iterativesolvers.html save_eigenspectrum_rdp private module subroutine save_eigenspectrum_rdp(lambda, residuals, fname)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: lambda(:)

Eigenalues.

real(kind=dp), intent(in) :: residuals(:)

Residual of the corresponding Ritz eigenpairs.

character(len=*), intent(in) :: fname

Name of the output file.

module procedure /github/workspace/API-doc/module/lightkrylov_iterativesolvers.html save_eigenspectrum_rsp private module subroutine save_eigenspectrum_rsp(lambda, residuals, fname)

Arguments

Type IntentOptional Attributes Name
real(kind=sp), intent(in) :: lambda(:)

Eigenalues.

real(kind=sp), intent(in) :: residuals(:)

Residual of the corresponding Ritz eigenpairs.

character(len=*), intent(in) :: fname

Name of the output file.