specialmatrices_symtridiagonal Module


Uses

    • stdlib_linalg_constants

Interfaces

public interface SymTridiagonal

This interface provides different methods to construct a SymTridiagonal matrix. Only the non-zero elements of are stored, i.e.

Syntax

  • Construct a SymTridiagonal matrix filled with zeros:
   integer, parameter :: n = 100
   type(SymTridiagonal) :: A

   A = SymTridiagonal(n)
  • Construct a SymTridiagonal matrix from rank-1 arrays:
   integer, parameter :: n
   real(dp), allocatable :: ev(:), dv(:)
   type(SymTridiagonal) :: A
   integer :: i

   dv = [(i, i=1, n)]; ev = [(2*i, i=1, n)]
   A = Tridiagonal(dv, ev)
  • Construct a SymTridiagonal matrix with constant diagonals:
   integer, parameter :: n
   real(dp), parameter :: d = 1.0_dp, e = 2.0_dp
   type(SymTridiagonal) :: A

   A = SymTridiagonal(d, e, n)

Note

Only double precision is currently supported for this matrix type.

Note

If is known to be symmetric positive definite, it can be constructed as A = SymTridiagonal(dv, ev, ifposdef=.true.):w

  • private pure module function construct(dv, ev, isposdef) result(A)

    Construct a SymTridiagonal matrix from the rank-1 arrays dv and ev.

    Arguments

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

    SymTridiagonal elements of the matrix.

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

    SymTridiagonal elements of the matrix.

    logical(kind=lk), intent(in), optional :: isposdef

    Whether A is positive-definite or not.

    Return Value type(SymTridiagonal)

    Symmetric Tridiagonal matrix.

  • private pure module function construct_constant(d, e, n, isposdef) result(A)

    Construct a SymTridiagonal matrix with constant diagonal elements.

    Arguments

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

    SymTridiagonal elements of the matrix.

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

    SymTridiagonal elements of the matrix.

    integer(kind=ilp), intent(in) :: n

    Dimension of the matrix.

    logical(kind=lk), intent(in), optional :: isposdef

    Whether A is positive-definite or not.

    Return Value type(SymTridiagonal)

    Symmetric Tridiagonal matrix.

  • private pure module function initialize(n) result(A)

    Construct a SymTridiagonal matrix filled with zeros.

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=ilp), intent(in) :: n

    Dimension of the matrix.

    Return Value type(SymTridiagonal)

    Symmetric Tridiagonal matrix.

public interface dense

This interface provides methods to convert a SymTridiagonal matrix to a regular rank-2 array.

Syntax

   B = dense(A)

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • B : Rank-2 array representation of the matrix .

  • private module function dense_rdp(A) result(B)

    Convert a SymTridiagonal matrix to a rank-2 array.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input diagonal matrix.

    Return Value real(kind=dp), allocatable, (:,:)

    Output dense rank-2 array.

public interface det

This interface overloads the det interface from stdlib_linag to compute the determinant where is of type SymTridiagonal.

Syntax

   d = det(A)

Arguments

  • A : Matrix of SymTridiagonal type. It is in an intent(in) argument.

  • d : Determinant of the matrix.

  • private pure module function det_rdp(A) result(d)

    Compute the determinant of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value real(kind=dp)

    Determinant of the matrix.

public interface eigh

This interface overloads the eigh interface from stdlib_linalg to compute the eigenvalues and eigenvectors of a matrix whose type is SymTridiagonal.

Syntax

   call eigh(A, lambda [, vectors])

Arguments

  • A : Matrix of SymTridiagonal. It is an intent(in) argument.

  • lambda : Rank-1 real array returning the eigenvalues of A in increasing order. It is an intent(out) argument.

  • vectors (optional) : Rank-2 array of the same kind as A returning the eigenvectors of A. It is an intent(out) argument.

  • private module subroutine eigh_rdp(A, lambda, vectors)

    Compute the eigenvalues and eigenvectors of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

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

    Eigenvalues.

    real(kind=dp), intent(out), optional, allocatable, target :: vectors(:,:)

    Eigenvectors.

public interface eigvalsh

This interface overloads the eigvalsh interface from stdlib_linalg to compute the eigenvalues of a matrix whose type is SymTridiagonal.

Syntax

   lambda = eigvalsh(A)

Arguments

  • A : real-valued matrix of SymTridiagonal type. It is an intent(in) argument.

  • lambda : Vector of eigenvalues in increasing order.

  • private module function eigvalsh_rdp(A) result(lambda)

    Utility function to compute the eigenvalues of a real SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value real(kind=dp), allocatable, (:)

    Eigenvalues.

public interface inv

  • private pure module function inv_rdp(A) result(B)

    Compute the inverse of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value real(kind=dp), allocatable, (:,:)

    Inverse of A.

public interface matmul

This interface overloads the Fortran intrinsic matmul for a SymTridiagonal matrix, both for matrix-vector and matrix-matrix products. For a matrix-matrix product , only the matrix has to be a SymTridiagonal matrix. Both and need to be standard Fortran rank-2 arrays. All the underlying functions are defined as pure.

Syntax

   y = matmul(A, x)
  • private module function spmv(A, x) result(y)

    Compute the matrix-vector product for a SymTridiagonal matrix . Both x and y are rank-1 arrays with the same kind as A.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    real(kind=dp), intent(in), target :: x(:)

    Input vector.

    Return Value real(kind=dp), target, allocatable, (:)

    Output vector.

  • private pure module function spmvs(A, X) result(Y)

    Compute the matrix-matrix product for a SymTridiagonal matrix and a dense matrix (rank-2 array). is also a rank-2 array with the same dimensions as .

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

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

    Input vectors.

    Return Value real(kind=dp), allocatable, (:,:)

    Output vectors.

public interface operator(*)

  • private pure module function scalar_multiplication_bis_rdp(A, alpha) result(B)

    Scalar multiplication with a SymTridiagonal matrix.

    Arguments

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

    Return Value type(SymTridiagonal)

  • private pure module function scalar_multiplication_rdp(alpha, A) result(B)

    Scalar multiplication with a SymTridiagonal matrix.

    Arguments

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

    Return Value type(SymTridiagonal)

public interface shape

  • private pure module function shape_rdp(A) result(arr_shape)

    Return the shape of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value integer(kind=ilp), (2)

    Shape of the matrix.

public interface size

  • private pure module function size_rdp(A, dim) result(arr_size)

    Return the size of SymTridiagonal matrix along a given dimension.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    integer(kind=ilp), intent(in), optional :: dim

    Queried dimension.

    Return Value integer(kind=ilp)

    Size of the matrix along the dimension dim.

public interface solve

This interface overloads the solve interface from stdlib_linalg for solving a linear system where is a SymTridiagonal matrix. It also enables to solve a linear system with multiple right-hand sides.

Syntax

   x = solve(A, b [, refine])

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • b : Rank-1 or rank-2 array defining the right-hand side(s). It is an intent(in) argument.

  • refine (optional) : Logical switch to enable solution refinement.

  • x : Solution of the linear system. It has the same type and shape as b.

  • private module function solve_multi_rhs(A, b, refine) result(x)

    Solve the linear system where is of type SymTridiagonal and B a standard rank-2 array. The solution matrix X has the same dimensions and kind as B.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Coefficient matrix.

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

    Right-hand side vectors.

    logical(kind=lk), intent(in), optional :: refine

    Whether iterative refinement of the solution is used or not.

    Return Value real(kind=dp), allocatable, (:,:)

    Solution vectors.

  • private module function solve_single_rhs(A, b, refine) result(x)

    Solve the linear system where is of type SymTridiagonal and b a standard rank-1 array. The solution vector x has the same dimension and kind as b.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Coefficient matrix.

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

    Right-hand side vector.

    logical(kind=lk), intent(in), optional :: refine

    Whether iterative refinement of the solution is used or not.

    Return Value real(kind=dp), allocatable, target, (:)

    Solution vector.

public interface svd

This interface overloads the svd interface from stdlib_linalg to compute the the singular value decomposition of a SymTridiagonal matrix .

Syntax

   call svd(A, s [, u] [, vt])

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • s : Rank-1 array real array returning the singular values of A. It is an intent(out) argument.

  • u (optional) : Rank-2 array of the same kind as A returning the left singular vectors of A as columns. Its size should be [n, n]. It is an intent(out) argument.

  • vt (optional) : Rank-2 array of the same kind as A returning the right singular vectors of A as rows. Its size should be [n, n]. It is an intent(out) argument.

  • private module subroutine svd_rdp(A, s, u, vt)

    Compute the singular value decomposition of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

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

    Singular values in descending order.

    real(kind=dp), intent(out), optional, allocatable :: u(:,:)

    Left singular vectors as columns.

    real(kind=dp), intent(out), optional, allocatable :: vt(:,:)

    Right singular vectors as rows.

public interface svdvals

This interface overloads the svdvals interface from stdlib_linalg to compute the singular values of a SymTridiagonal matrix .

Syntax

   s = svdvals(A)

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • s : Vector of singular values sorted in decreasing order.

  • private module function svdvals_rdp(A) result(s)

    Compute the singular values of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value real(kind=dp), allocatable, (:)

    Singular values in descending order.

public interface trace

This interface overloads the trace interface from stdlib_linalg to compute the trace of a matrix of type SymTridiagonal.

Syntax

   tr = trace(A)

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • tr: Trace of the matrix.

  • private pure module function trace_rdp(A) result(tr)

    Compute the trace of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value real(kind=dp)

    Trace of the matrix.

public interface transpose

This interface overloads the Fortran intrinsic procedure to define the transpose operation for a SymTridiagonal matrix.

Syntax

   B = transpose(A)

Arguments

  • A : Matrix of SymTridiagonal type. It is an intent(in) argument.

  • B : Resulting transposed matrix. It is of the same type as A.

  • private pure module function transpose_rdp(A) result(B)

    Compute the transpose of a SymTridiagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    type(SymTridiagonal), intent(in) :: A

    Input matrix.

    Return Value type(SymTridiagonal)

    Transpose of the matrix.


Derived Types

type, public ::  SymTridiagonal

Base type used to define a SymTridiagonal matrix of size [n, n] with diagonals given by rank-1 arrays dv (size n) and ev (size n-1).

Constructor

This interface provides different methods to construct a SymTridiagonal matrix. Only the non-zero elements of are stored, i.e.

Read more…
private pure, module function construct (dv, ev, isposdef)

Construct a SymTridiagonal matrix from the rank-1 arrays dv and ev.

private pure, module function construct_constant (d, e, n, isposdef)

Construct a SymTridiagonal matrix with constant diagonal elements.

private pure, module function initialize (n)

Construct a SymTridiagonal matrix filled with zeros.