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

Eigendecomposition of a square matrix, returning eigenvalues and optionally eigenvectors. More...

Public Member Functions

subroutine la_eig_standard_s (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_standard_d (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_standard_q (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_standard_c (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_standard_z (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_standard_w (a, lambda, right, left, overwrite_a, err)
 
subroutine la_real_eig_standard_s (a, lambda, right, left, overwrite_a, err)
 
subroutine la_real_eig_standard_d (a, lambda, right, left, overwrite_a, err)
 
subroutine la_real_eig_standard_q (a, lambda, right, left, overwrite_a, err)
 
subroutine la_eig_generalized_s (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_eig_generalized_d (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_eig_generalized_q (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_eig_generalized_c (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_eig_generalized_z (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_eig_generalized_w (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_real_eig_generalized_s (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_real_eig_generalized_d (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 
subroutine la_real_eig_generalized_q (a, b, lambda, right, left, overwrite_a, overwrite_b, err)
 

Detailed Description

Eigendecomposition of a square matrix, returning eigenvalues and optionally eigenvectors.

Summary

This interface provides methods for computing the eigenvalues and eigenvectors of a real or complex matrix. It supports both standard and generalized eigenvalue problems, allowing for the decomposition of a matrix A alone or a pair of matrices (A, B) in the generalized case.

Description

Given a square matrix A , this routine computes its eigenvalues \lambda and, optionally, its right or left eigenvectors:

A v = \lambda v,

where v represents an eigenvector corresponding to eigenvalue \lambda .

In the generalized eigenvalue problem case, the routine solves:

[ A v = \lambda B v. ]

The computation supports both real and complex matrices. If requested, eigenvectors are returned as additional output arguments. The function provides options to allow in-place modification of A and B for performance optimization.

Note
The solution is based on LAPACK's GEEV and GGEV routines.

Arguments

  • a: A real or complex matrix of size [n,n] , representing the input matrix to be decomposed.
  • b (optional): A matrix of size [n,n] , same type and kind as a, representing the second matrix in the generalized eigenvalue problem.
  • lambda: A complex or real array of length n containing the computed eigenvalues.
  • right (optional): A complex matrix of size [n,n] containing the right eigenvectors as columns.
  • left (optional): A complex matrix of size [n,n] containing the left eigenvectors as columns.
  • overwrite_a (optional): A logical flag indicating whether A can be overwritten for performance optimization.
  • overwrite_b (optional): A logical flag indicating whether B can be overwritten (only in the generalized case).
  • err (optional): A la_state variable to handle errors. If not provided, execution will stop on errors.

Errors

  • Raises LINALG_VALUE_ERROR if the matrices have invalid/incompatible sizes.
  • Raises LINALG_ERROR if the eigendecomposition fails.
  • If err is not provided, exceptions will trigger an error stop.

Notes

  • The computed eigenvectors are normalized.
  • If computing real eigenvalues, an error is returned if eigenvalues have nonzero imaginary parts.
  • This routine is based on LAPACK's GEEV and GGEV routines.
  • Overwriting A or B can improve performance but destroys the original matrix data.

Member Function/Subroutine Documentation

◆ la_eig_generalized_c()

subroutine la_eig::eig::la_eig_generalized_c ( complex(sp), dimension(:,:), intent(inout), target a,
complex(sp), dimension(:,:), intent(inout), target b,
complex(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_generalized_d()

subroutine la_eig::eig::la_eig_generalized_d ( real(dp), dimension(:,:), intent(inout), target a,
real(dp), dimension(:,:), intent(inout), target b,
complex(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_generalized_q()

subroutine la_eig::eig::la_eig_generalized_q ( real(qp), dimension(:,:), intent(inout), target a,
real(qp), dimension(:,:), intent(inout), target b,
complex(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_generalized_s()

subroutine la_eig::eig::la_eig_generalized_s ( real(sp), dimension(:,:), intent(inout), target a,
real(sp), dimension(:,:), intent(inout), target b,
complex(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_generalized_w()

subroutine la_eig::eig::la_eig_generalized_w ( complex(qp), dimension(:,:), intent(inout), target a,
complex(qp), dimension(:,:), intent(inout), target b,
complex(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_generalized_z()

subroutine la_eig::eig::la_eig_generalized_z ( complex(dp), dimension(:,:), intent(inout), target a,
complex(dp), dimension(:,:), intent(inout), target b,
complex(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_c()

subroutine la_eig::eig::la_eig_standard_c ( complex(sp), dimension(:,:), intent(inout), target a,
complex(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_d()

subroutine la_eig::eig::la_eig_standard_d ( real(dp), dimension(:,:), intent(inout), target a,
complex(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_q()

subroutine la_eig::eig::la_eig_standard_q ( real(qp), dimension(:,:), intent(inout), target a,
complex(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_s()

subroutine la_eig::eig::la_eig_standard_s ( real(sp), dimension(:,:), intent(inout), target a,
complex(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_w()

subroutine la_eig::eig::la_eig_standard_w ( complex(qp), dimension(:,:), intent(inout), target a,
complex(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_eig_standard_z()

subroutine la_eig::eig::la_eig_standard_z ( complex(dp), dimension(:,:), intent(inout), target a,
complex(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of eigenvalues
[out]right[optional] RIGHT eigenvectors of A (as columns)
[out]left[optional] LEFT eigenvectors of A (as columns)
[in]overwrite_a[optional] Can A data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_generalized_d()

subroutine la_eig::eig::la_real_eig_generalized_d ( real(dp), dimension(:,:), intent(inout), target a,
real(dp), dimension(:,:), intent(inout), target b,
real(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_generalized_q()

subroutine la_eig::eig::la_real_eig_generalized_q ( real(qp), dimension(:,:), intent(inout), target a,
real(qp), dimension(:,:), intent(inout), target b,
real(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_generalized_s()

subroutine la_eig::eig::la_real_eig_generalized_s ( real(sp), dimension(:,:), intent(inout), target a,
real(sp), dimension(:,:), intent(inout), target b,
real(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional overwrite_b,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[in,out]bGeneralized problem matrix B[n,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]overwrite_b[optional] Can B data be overwritten and destroyed? (default: no)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_standard_d()

subroutine la_eig::eig::la_real_eig_standard_d ( real(dp), dimension(:,:), intent(inout), target a,
real(dp), dimension(:), intent(out) lambda,
complex(dp), dimension(:,:), intent(out), optional, target right,
complex(dp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_standard_q()

subroutine la_eig::eig::la_real_eig_standard_q ( real(qp), dimension(:,:), intent(inout), target a,
real(qp), dimension(:), intent(out) lambda,
complex(qp), dimension(:,:), intent(out), optional, target right,
complex(qp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_real_eig_standard_s()

subroutine la_eig::eig::la_real_eig_standard_s ( real(sp), dimension(:,:), intent(inout), target a,
real(sp), dimension(:), intent(out) lambda,
complex(sp), dimension(:,:), intent(out), optional, target right,
complex(sp), dimension(:,:), intent(out), optional, target left,
logical(lk), intent(in), optional overwrite_a,
type(la_state), intent(out), optional err )
Parameters
[in,out]aInput matrix A[m,n]
[out]lambdaArray of real eigenvalues
[out]rightThe columns of RIGHT contain the right eigenvectors of A
[out]leftThe columns of LEFT contain the left eigenvectors of A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[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: