Skip to content

Commit e6aa0f5

Browse files
committed
Option to output real eigenvalues only
1 parent 36c8d5a commit e6aa0f5

File tree

3 files changed

+64
-2
lines changed

3 files changed

+64
-2
lines changed

Diff for: doc/specs/stdlib_linalg.md

+2-1
Original file line numberDiff line numberDiff line change
@@ -923,7 +923,7 @@ The solver is based on LAPACK's `*GEEV` backends.
923923

924924
`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.
925925

926-
`lambda`: Shall be a `complex` rank-1 array of the same kind as `a`, containing the eigenvalues. It is an `intent(out)` argument.
926+
`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.
927927

928928
`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.
929929

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

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

942943
### Example

Diff for: src/stdlib_linalg.fypp

+21
Original file line numberDiff line numberDiff line change
@@ -597,6 +597,27 @@ module stdlib_linalg
597597
end subroutine stdlib_linalg_eig_${ri}$
598598
#:endif
599599
#:endfor
600+
#:for rk,rt,ri in REAL_KINDS_TYPES
601+
#:if rk!="xdp"
602+
module subroutine stdlib_linalg_real_eig_${ri}$(a,lambda,right,left,overwrite_a,err)
603+
!! Eigendecomposition of matrix A returning an array `lambda` of real eigenvalues,
604+
!! and optionally right or left eigenvectors. Returns an error if the eigenvalues had
605+
!! non-trivial imaginary parts.
606+
!> Input matrix A[m,n]
607+
${rt}$, intent(inout), target :: a(:,:)
608+
!> Array of real eigenvalues
609+
real(${rk}$), intent(out) :: lambda(:)
610+
!> The columns of RIGHT contain the right eigenvectors of A
611+
complex(${rk}$), optional, intent(out), target :: right(:,:)
612+
!> The columns of LEFT contain the left eigenvectors of A
613+
complex(${rk}$), optional, intent(out), target :: left(:,:)
614+
!> [optional] Can A data be overwritten and destroyed?
615+
logical(lk), optional, intent(in) :: overwrite_a
616+
!> [optional] state return flag. On error if not requested, the code will stop
617+
type(linalg_state_type), optional, intent(out) :: err
618+
end subroutine stdlib_linalg_real_eig_${ri}$
619+
#:endif
620+
#:endfor
600621
end interface eig
601622

602623
! Eigenvalues of a square matrix

Diff for: src/stdlib_linalg_eigenvalues.fypp

+41-1
Original file line numberDiff line numberDiff line change
@@ -315,6 +315,47 @@ submodule (stdlib_linalg) stdlib_linalg_eigenvalues
315315

316316
end subroutine stdlib_linalg_eig_${ri}$
317317

318+
module subroutine stdlib_linalg_real_eig_${ri}$(a,lambda,right,left,overwrite_a,err)
319+
!! Eigendecomposition of matrix A returning an array `lambda` of real eigenvalues,
320+
!! and optionally right or left eigenvectors. Returns an error if the eigenvalues had
321+
!! non-trivial imaginary parts.
322+
!> Input matrix A[m,n]
323+
${rt}$, intent(inout), target :: a(:,:)
324+
!> Array of real eigenvalues
325+
real(${rk}$), intent(out) :: lambda(:)
326+
!> The columns of RIGHT contain the right eigenvectors of A
327+
complex(${rk}$), optional, intent(out), target :: right(:,:)
328+
!> The columns of LEFT contain the left eigenvectors of A
329+
complex(${rk}$), optional, intent(out), target :: left(:,:)
330+
!> [optional] Can A data be overwritten and destroyed?
331+
logical(lk), optional, intent(in) :: overwrite_a
332+
!> [optional] state return flag. On error if not requested, the code will stop
333+
type(linalg_state_type), optional, intent(out) :: err
334+
335+
type(linalg_state_type) :: err0
336+
integer(ilp) :: n
337+
complex(${rk}$), allocatable :: clambda(:)
338+
real(${rk}$), parameter :: rtol = epsilon(0.0_${rk}$)
339+
real(${rk}$), parameter :: atol = tiny(0.0_${rk}$)
340+
341+
n = size(lambda,dim=1,kind=ilp)
342+
allocate(clambda(n))
343+
344+
call stdlib_linalg_eig_${ri}$(a,clambda,right,left,overwrite_a,err0)
345+
346+
! Check that no eigenvalues have meaningful imaginary part
347+
if (err0%ok() .and. any(aimag(clambda)>atol+rtol*abs(abs(clambda)))) then
348+
err0 = linalg_state_type(this,LINALG_VALUE_ERROR, &
349+
'complex eigenvalues detected: max(imag(lambda))=',maxval(aimag(clambda)))
350+
endif
351+
352+
! Return real components only
353+
lambda(:n) = real(clambda,kind=${rk}$)
354+
355+
call linalg_error_handling(err0,err)
356+
357+
end subroutine stdlib_linalg_real_eig_${ri}$
358+
318359
module function stdlib_linalg_eigvalsh_${ri}$(a,upper_a,err) result(lambda)
319360
!! Return an array of eigenvalues of real symmetric / complex hermitian A
320361
!> Input matrix A[m,n]
@@ -535,5 +576,4 @@ submodule (stdlib_linalg) stdlib_linalg_eigenvalues
535576
#:endif
536577
#:endfor
537578

538-
539579
end submodule stdlib_linalg_eigenvalues

0 commit comments

Comments
 (0)