Skip to content

Commit

Permalink
Option to output real eigenvalues only
Browse files Browse the repository at this point in the history
  • Loading branch information
perazz committed Jun 5, 2024
1 parent 36c8d5a commit e6aa0f5
Show file tree
Hide file tree
Showing 3 changed files with 64 additions and 2 deletions.
3 changes: 2 additions & 1 deletion doc/specs/stdlib_linalg.md
Original file line number Diff line number Diff line change
Expand Up @@ -923,7 +923,7 @@ The solver is based on LAPACK's `*GEEV` backends.

`a` : `real` or `complex` square array containing the coefficient matrix. If `overwrite_a=.false.`, it is an `intent(in)` argument. Otherwise, it is an `intent(inout)` argument and is destroyed by the call.

`lambda`: Shall be a `complex` rank-1 array of the same kind as `a`, containing the eigenvalues. It is an `intent(out)` argument.
`lambda`: Shall be a `complex` or `real` rank-1 array of the same kind as `a`, containing the eigenvalues, or their `real` component only. It is an `intent(out)` argument.

`right` (optional): Shall be a `complex` rank-2 array of the same size and kind as `a`, containing the right eigenvectors of `a`. It is an `intent(out)` argument.

Expand All @@ -937,6 +937,7 @@ The solver is based on LAPACK's `*GEEV` backends.

Raises `LINALG_ERROR` if the calculation did not converge.
Raises `LINALG_VALUE_ERROR` if any matrix or arrays have invalid/incompatible sizes.
Raises `LINALG_VALUE_ERROR` if the `real` component is only requested, but the eigenvalues have non-trivial imaginary parts.
If `err` is not present, exceptions trigger an `error stop`.

### Example
Expand Down
21 changes: 21 additions & 0 deletions src/stdlib_linalg.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -597,6 +597,27 @@ module stdlib_linalg
end subroutine stdlib_linalg_eig_${ri}$
#:endif
#:endfor
#:for rk,rt,ri in REAL_KINDS_TYPES
#:if rk!="xdp"
module subroutine stdlib_linalg_real_eig_${ri}$(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.
!> Input matrix A[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> Array of real eigenvalues
real(${rk}$), intent(out) :: lambda(:)
!> The columns of RIGHT contain the right eigenvectors of A
complex(${rk}$), optional, intent(out), target :: right(:,:)
!> The columns of LEFT contain the left eigenvectors of A
complex(${rk}$), optional, intent(out), target :: left(:,:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err
end subroutine stdlib_linalg_real_eig_${ri}$
#:endif
#:endfor
end interface eig
! Eigenvalues of a square matrix
Expand Down
42 changes: 41 additions & 1 deletion src/stdlib_linalg_eigenvalues.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -315,6 +315,47 @@ submodule (stdlib_linalg) stdlib_linalg_eigenvalues

end subroutine stdlib_linalg_eig_${ri}$

module subroutine stdlib_linalg_real_eig_${ri}$(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.
!> Input matrix A[m,n]
${rt}$, intent(inout), target :: a(:,:)
!> Array of real eigenvalues
real(${rk}$), intent(out) :: lambda(:)
!> The columns of RIGHT contain the right eigenvectors of A
complex(${rk}$), optional, intent(out), target :: right(:,:)
!> The columns of LEFT contain the left eigenvectors of A
complex(${rk}$), optional, intent(out), target :: left(:,:)
!> [optional] Can A data be overwritten and destroyed?
logical(lk), optional, intent(in) :: overwrite_a
!> [optional] state return flag. On error if not requested, the code will stop
type(linalg_state_type), optional, intent(out) :: err

type(linalg_state_type) :: err0
integer(ilp) :: n
complex(${rk}$), allocatable :: clambda(:)
real(${rk}$), parameter :: rtol = epsilon(0.0_${rk}$)
real(${rk}$), parameter :: atol = tiny(0.0_${rk}$)

n = size(lambda,dim=1,kind=ilp)
allocate(clambda(n))

call stdlib_linalg_eig_${ri}$(a,clambda,right,left,overwrite_a,err0)

! Check that no eigenvalues have meaningful imaginary part
if (err0%ok() .and. any(aimag(clambda)>atol+rtol*abs(abs(clambda)))) then
err0 = linalg_state_type(this,LINALG_VALUE_ERROR, &
'complex eigenvalues detected: max(imag(lambda))=',maxval(aimag(clambda)))
endif

! Return real components only
lambda(:n) = real(clambda,kind=${rk}$)

call linalg_error_handling(err0,err)

end subroutine stdlib_linalg_real_eig_${ri}$

module function stdlib_linalg_eigvalsh_${ri}$(a,upper_a,err) result(lambda)
!! Return an array of eigenvalues of real symmetric / complex hermitian A
!> Input matrix A[m,n]
Expand Down Expand Up @@ -535,5 +576,4 @@ submodule (stdlib_linalg) stdlib_linalg_eigenvalues
#:endif
#:endfor


end submodule stdlib_linalg_eigenvalues

0 comments on commit e6aa0f5

Please sign in to comment.