Solves the eigendecomposition for a real symmetric or complex Hermitian square matrix. (Specification)
Subroutine interface for computing eigenvalues and eigenvectors of a real symmetric or complex Hermitian square matrix.
This interface provides methods for computing the eigenvalues, and optionally eigenvectors,
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). The vectors are orthogonal, and may be returned as columns of an optional
matrix vectors
with the same kind and size as A
.
Preallocated space for both eigenvalues lambda
and the eigenvector matrix must be user-provided.
Note
The solution is based on LAPACK's eigenproblem solvers *SYEV
/*HEEV
.
Note
BLAS/LAPACK backends do not currently support extended precision (xdp
).
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda
of eigenvalues, and optionally right or left eigenvectors.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=sp), | intent(inout), | target | :: | a(:,:) |
Input matrix A[m,n] |
|
real(kind=sp), | intent(out) | :: | lambda(:) |
Array of eigenvalues |
||
complex(kind=sp), | intent(out), | optional, | target | :: | vectors(:,:) |
The columns of vectors contain the orthonormal eigenvectors of A |
logical(kind=lk), | intent(in), | optional | :: | upper_a |
[optional] Should the upper/lower half of A be used? Default: lower |
|
logical(kind=lk), | intent(in), | optional | :: | overwrite_a |
[optional] Can A data be overwritten and destroyed? |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda
of eigenvalues, and optionally right or left eigenvectors.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(inout), | target | :: | a(:,:) |
Input matrix A[m,n] |
|
real(kind=dp), | intent(out) | :: | lambda(:) |
Array of eigenvalues |
||
real(kind=dp), | intent(out), | optional, | target | :: | vectors(:,:) |
The columns of vectors contain the orthonormal eigenvectors of A |
logical(kind=lk), | intent(in), | optional | :: | upper_a |
[optional] Should the upper/lower half of A be used? Default: lower |
|
logical(kind=lk), | intent(in), | optional | :: | overwrite_a |
[optional] Can A data be overwritten and destroyed? |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda
of eigenvalues, and optionally right or left eigenvectors.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=sp), | intent(inout), | target | :: | a(:,:) |
Input matrix A[m,n] |
|
real(kind=sp), | intent(out) | :: | lambda(:) |
Array of eigenvalues |
||
real(kind=sp), | intent(out), | optional, | target | :: | vectors(:,:) |
The columns of vectors contain the orthonormal eigenvectors of A |
logical(kind=lk), | intent(in), | optional | :: | upper_a |
[optional] Should the upper/lower half of A be used? Default: lower |
|
logical(kind=lk), | intent(in), | optional | :: | overwrite_a |
[optional] Can A data be overwritten and destroyed? |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda
of eigenvalues, and optionally right or left eigenvectors.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=dp), | intent(inout), | target | :: | a(:,:) |
Input matrix A[m,n] |
|
real(kind=dp), | intent(out) | :: | lambda(:) |
Array of eigenvalues |
||
complex(kind=dp), | intent(out), | optional, | target | :: | vectors(:,:) |
The columns of vectors contain the orthonormal eigenvectors of A |
logical(kind=lk), | intent(in), | optional | :: | upper_a |
[optional] Should the upper/lower half of A be used? Default: lower |
|
logical(kind=lk), | intent(in), | optional | :: | overwrite_a |
[optional] Can A data be overwritten and destroyed? |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |