Processing math: 100%
fortran-lapack
All Classes Namespaces Files Functions Variables Pages
la_eig::eigh Interface Reference

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.
 

Detailed Description

Eigendecomposition of a real symmetric or complex Hermitian matrix.

Summary

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.

Description

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.

Note
This method is based on LAPACK's SYEVD routine for real symmetric matrices and HEEVD routine for complex Hermitian matrices.

Arguments

  • 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.

Return value

  • lambda: A real array containing the eigenvalues of matrix A , computed as part of the decomposition.

Errors

  • Raises LINALG_VALUE_ERROR if the matrix is not symmetric or Hermitian.
  • Raises LINALG_ERROR if the eigendecomposition fails.
  • If err is not provided, execution will stop on errors.

Notes

  • The eigenvalues are always real for symmetric and Hermitian matrices.
  • If the matrix is Hermitian, the result will be complex eigenvectors even if the input matrix is real.
  • This method is based on LAPACK's SYEVD for real symmetric matrices and HEEVD for complex Hermitian matrices.
  • The matrix A can be overwritten in-place to improve performance but its original data will be lost.

Member Function/Subroutine Documentation

◆ la_eigh_c()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

◆ la_eigh_d()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

◆ la_eigh_q()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

◆ la_eigh_s()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

◆ la_eigh_w()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

◆ la_eigh_z()

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.

Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]vectorsThe 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

The documentation for this interface was generated from the following file: