eig Interface

public interface eig

Solves the eigendecomposition for square matrix . (Specification)

Summary

Subroutine interface for computing eigenvalues and eigenvectors of a square matrix.

Description

This interface provides methods for computing the eigenvalues, and optionally eigenvectors, of a general square matrix. Supported data types include real and complex, and no assumption is made on the matrix structure. The user may request either left, right, or both eigenvectors to be returned. They are returned as columns of a square matrix with the same size as A. Preallocated space for both eigenvalues lambda and the eigenvector matrices must be user-provided.

Note

The solution is based on LAPACK's general eigenproblem solvers *GEEV.

Note

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


Subroutines

private module subroutine stdlib_linalg_eig_c(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.

Arguments

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

Input matrix A[m,n]

complex(kind=sp), intent(out) :: lambda(:)

Array of eigenvalues

complex(kind=sp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=sp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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

private module subroutine stdlib_linalg_eig_d(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.

Arguments

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

Input matrix A[m,n]

complex(kind=dp), intent(out) :: lambda(:)

Array of eigenvalues

complex(kind=dp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=dp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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

private module subroutine stdlib_linalg_eig_s(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.

Arguments

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

Input matrix A[m,n]

complex(kind=sp), intent(out) :: lambda(:)

Array of eigenvalues

complex(kind=sp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=sp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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

private module subroutine stdlib_linalg_eig_z(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of eigenvalues, and optionally right or left eigenvectors.

Arguments

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

Input matrix A[m,n]

complex(kind=dp), intent(out) :: lambda(:)

Array of eigenvalues

complex(kind=dp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=dp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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

private module subroutine stdlib_linalg_real_eig_d(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of real eigenvalues, and optionally right or left eigenvectors. Returns an error if the eigenvalues had non-trivial imaginary parts.

Arguments

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

Input matrix A[m,n]

real(kind=dp), intent(out) :: lambda(:)

Array of real eigenvalues

complex(kind=dp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=dp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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

private module subroutine stdlib_linalg_real_eig_s(a, lambda, right, left, overwrite_a, err)

Eigendecomposition of matrix A returning an array lambda of real eigenvalues, and optionally right or left eigenvectors. Returns an error if the eigenvalues had non-trivial imaginary parts.

Arguments

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

Input matrix A[m,n]

real(kind=sp), intent(out) :: lambda(:)

Array of real eigenvalues

complex(kind=sp), intent(out), optional, target :: right(:,:)

The columns of RIGHT contain the right eigenvectors of A

complex(kind=sp), intent(out), optional, target :: left(:,:)

The columns of LEFT contain the left eigenvectors of A

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