eigvalsh Interface

public interface eigvalsh

Returns the eigenvalues , , for a real symmetric or complex Hermitian square matrix. (Specification)

Summary

Function interface for computing the eigenvalues of a real symmetric or complex hermitian square matrix.

Description

This interface provides functions for returning the eigenvalues of a real symmetric or complex Hermitian square matrix. Supported data types include real and complex. The matrix must be symmetric (if real) or Hermitian (if complex). Only the lower or upper half of the matrix is accessed, and the user can select which using the optional upper_a flag (default: use lower half). An error stop is thrown in case of failure; otherwise, error information can be returned as an optional type(linalg_state_type) output flag.

Note

The solution is based on LAPACK's eigenproblem solvers *SYEV/*HEEV.

Note

BLAS/LAPACK backends do not currently support extended precision (xdp).


Functions

private module function stdlib_linalg_eigvalsh_c(a, upper_a, err) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

type(linalg_state_type), intent(out) :: err

[optional] state return flag. On error if not requested, the code will stop

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

Array of singular values

private module function stdlib_linalg_eigvalsh_d(a, upper_a, err) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

type(linalg_state_type), intent(out) :: err

[optional] state return flag. On error if not requested, the code will stop

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

Array of singular values

private module function stdlib_linalg_eigvalsh_noerr_c(a, upper_a) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

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

Array of singular values

private module function stdlib_linalg_eigvalsh_noerr_d(a, upper_a) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

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

Array of singular values

private module function stdlib_linalg_eigvalsh_noerr_s(a, upper_a) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

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

Array of singular values

private module function stdlib_linalg_eigvalsh_noerr_z(a, upper_a) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

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

Array of singular values

private module function stdlib_linalg_eigvalsh_s(a, upper_a, err) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

type(linalg_state_type), intent(out) :: err

[optional] state return flag. On error if not requested, the code will stop

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

Array of singular values

private module function stdlib_linalg_eigvalsh_z(a, upper_a, err) result(lambda)

Return an array of eigenvalues of real symmetric / complex hermitian A

Arguments

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

Input matrix A[m,n]

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

[optional] Should the upper/lower half of A be used? Default: lower

type(linalg_state_type), intent(out) :: err

[optional] state return flag. On error if not requested, the code will stop

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

Array of singular values