|
fortran-lapack
|
Eigendecomposition of a real symmetric or complex Hermitian matrix. More...
Public Member Functions | |
| subroutine | la_eigh_s (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
| subroutine | la_eigh_d (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
| subroutine | la_eigh_q (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
| subroutine | la_eigh_c (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
| subroutine | la_eigh_z (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
| subroutine | la_eigh_w (a, lambda, vectors, upper_a, overwrite_a, err) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors. | |
Eigendecomposition of a real symmetric or complex Hermitian matrix.
This interface provides methods for computing the eigenvalues and optionally the eigenvectors of a real symmetric or complex Hermitian matrix. The decomposition computes the eigenvalues of matrix \( A \) and, optionally, the eigenvectors corresponding to those eigenvalues. The routine is designed to handle both real and complex matrices.
Given a real symmetric or complex Hermitian matrix \( A \), this routine computes its eigenvalues \( \lambda \) and, optionally, the right or left eigenvectors:
\[ A v = \lambda v \]
where \( v \) represents an eigenvector corresponding to eigenvalue \( \lambda \). The eigenvalues are always returned as real numbers, and the eigenvectors, if requested, are orthonormal.
This function uses the LAPACK routines designed for symmetric or Hermitian matrices, ensuring numerical stability and performance.
By default, the routine uses the lower half of the matrix for the computation. However, the upper half can optionally be used by setting the upper_a argument.
a: A real or complex symmetric or Hermitian matrix of size \( [n,n] \), representing the matrix to be decomposed. On input, it contains the matrix \( A \), and on output, it is overwritten with the eigenvalues in the diagonal.lambda: A real array of length \( n \), containing the computed eigenvalues of \( A \).vectors (optional): A real or complex matrix of size \( [n,n] \) containing the orthonormal eigenvectors of \( A \), with the eigenvectors stored as columns.overwrite_a (optional): A logical flag indicating whether the matrix \( A \) can be overwritten for performance optimization.upper_a (optional): A logical flag indicating whether the upper half of \( A \) should be used (default is the lower half).err (optional): A la_state variable to handle errors. If not provided, execution will stop on errors.lambda: A real array containing the eigenvalues of matrix \( A \), computed as part of the decomposition.err is not provided, execution will stop on errors.| subroutine la_eig::eigh::la_eigh_c | ( | complex(sp), dimension(:,:), intent(inout), target | a, |
| real(sp), dimension(:), intent(out) | lambda, | ||
| complex(sp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |
| subroutine la_eig::eigh::la_eigh_d | ( | real(dp), dimension(:,:), intent(inout), target | a, |
| real(dp), dimension(:), intent(out) | lambda, | ||
| real(dp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |
| subroutine la_eig::eigh::la_eigh_q | ( | real(qp), dimension(:,:), intent(inout), target | a, |
| real(qp), dimension(:), intent(out) | lambda, | ||
| real(qp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |
| subroutine la_eig::eigh::la_eigh_s | ( | real(sp), dimension(:,:), intent(inout), target | a, |
| real(sp), dimension(:), intent(out) | lambda, | ||
| real(sp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |
| subroutine la_eig::eigh::la_eigh_w | ( | complex(qp), dimension(:,:), intent(inout), target | a, |
| real(qp), dimension(:), intent(out) | lambda, | ||
| complex(qp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |
| subroutine la_eig::eigh::la_eigh_z | ( | complex(dp), dimension(:,:), intent(inout), target | a, |
| real(dp), dimension(:), intent(out) | lambda, | ||
| complex(dp), dimension(:,:), intent(out), optional, target | vectors, | ||
| logical(lk), intent(in), optional | upper_a, | ||
| logical(lk), intent(in), optional | overwrite_a, | ||
| type(la_state), intent(out), optional | err ) |
Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.
| [in,out] | a | Input matrix A[m,n] |
| [out] | lambda | Array of eigenvalues |
| [out] | vectors | The columns of vectors contain the orthonormal eigenvectors of A |
| [in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
| [in] | upper_a | [optional] Should the upper/lower half of A be used? Default: lower |
| [out] | err | [optional] state return flag. On error if not requested, the code will stop |