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 by the AbstractLinops module. It is an intent(in) 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 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.

Note

Although the interface to pass a preconditioner is exposed, it is not currently implemented.

  • public subroutine cg_rsp(A, b, x, info, rtol, atol, preconditioner, options)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rsp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (not yet supported).

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

    Options for the conjugate gradient solver.

  • public subroutine cg_rdp(A, b, x, info, rtol, atol, preconditioner, options)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rdp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (not yet supported).

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

    Options for the conjugate gradient solver.

  • public subroutine cg_csp(A, b, x, info, rtol, atol, preconditioner, options)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_csp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (not yet supported).

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

    Options for the conjugate gradient solver.

  • public subroutine cg_cdp(A, b, x, info, rtol, atol, preconditioner, options)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_cdp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (not yet supported).

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

    Options for the conjugate gradient solver.

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(in) 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.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rsp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_sym_linop_rdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_csp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_hermitian_linop_cdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance

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(in) 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, kdim, select, tolerance, transpose)

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rsp), intent(in) :: 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.

    integer, intent(in), optional :: kdim
    procedure(eigvals_select_sp), optional :: select

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim
    procedure(eigvals_select_dp), optional :: select

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_csp), intent(in) :: 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.

    integer, intent(in), optional :: kdim
    procedure(eigvals_select_sp), optional :: select

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_cdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim
    procedure(eigvals_select_dp), optional :: select

    Desired number of eigenpairs.

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

    Tolerance.

    logical, intent(in), optional :: transpose

    Determine whether or is being used.

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(in) 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.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rsp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (optional).

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

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rdp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (optional).

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

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_csp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (optional).

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

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_cdp), intent(in) :: 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(in), optional :: preconditioner

    Preconditioner (optional).

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

    GMRES options.

    logical, intent(in), optional :: transpose

    Whether or is being used.

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 subroutine save_eigenspectrum_sp(eigvals, 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) :: eigvals(:)

    Eigenalues.

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

    Residual of the corresponding Ritz eigenpairs.

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

    Name of the output file.

  • private subroutine save_eigenspectrum_dp(eigvals, 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) :: eigvals(:)

    Eigenalues.

    real(kind=dp), 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(in) 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_sp or tolerance = rtol_dp.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rsp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_csp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_cdp), intent(in) :: 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.

    integer, intent(in), optional :: kdim

    Desired number of eigenpairs.

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

    Tolerance.


Abstract Interfaces

abstract interface

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

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_cdp), intent(in) :: 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(in), 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.

abstract interface

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

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_csp), intent(in) :: 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(in), 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.

abstract interface

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

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rdp), intent(in) :: 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(in), 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.

abstract interface

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

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

    Arguments

    Type IntentOptional Attributes Name
    class(abstract_linop_rsp), intent(in) :: 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(in), 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.


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

Subroutines

public subroutine cg_cdp(A, b, x, info, rtol, atol, preconditioner, options)

Arguments

Type IntentOptional Attributes Name
class(abstract_hermitian_linop_cdp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (not yet supported).

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

Options for the conjugate gradient solver.

public subroutine cg_csp(A, b, x, info, rtol, atol, preconditioner, options)

Arguments

Type IntentOptional Attributes Name
class(abstract_hermitian_linop_csp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (not yet supported).

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

Options for the conjugate gradient solver.

public subroutine cg_rdp(A, b, x, info, rtol, atol, preconditioner, options)

Arguments

Type IntentOptional Attributes Name
class(abstract_sym_linop_rdp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (not yet supported).

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

Options for the conjugate gradient solver.

public subroutine cg_rsp(A, b, x, info, rtol, atol, preconditioner, options)

Arguments

Type IntentOptional Attributes Name
class(abstract_sym_linop_rsp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (not yet supported).

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

Options for the conjugate gradient solver.

public subroutine gmres_cdp(A, b, x, info, rtol, atol, preconditioner, options, transpose)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_cdp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (optional).

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

GMRES options.

logical, intent(in), optional :: transpose

Whether or is being used.

public subroutine gmres_csp(A, b, x, info, rtol, atol, preconditioner, options, transpose)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_csp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (optional).

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

GMRES options.

logical, intent(in), optional :: transpose

Whether or is being used.

public subroutine gmres_rdp(A, b, x, info, rtol, atol, preconditioner, options, transpose)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_rdp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (optional).

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

GMRES options.

logical, intent(in), optional :: transpose

Whether or is being used.

public subroutine gmres_rsp(A, b, x, info, rtol, atol, preconditioner, options, transpose)

Arguments

Type IntentOptional Attributes Name
class(abstract_linop_rsp), intent(in) :: 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(in), optional :: preconditioner

Preconditioner (optional).

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

GMRES options.

logical, intent(in), optional :: transpose

Whether or is being used.