Skip to content

Commit 6ce5172

Browse files
committed
submodule
1 parent 7645a45 commit 6ce5172

File tree

2 files changed

+112
-64
lines changed

2 files changed

+112
-64
lines changed

src/stdlib_linalg.fypp

Lines changed: 102 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -552,6 +552,108 @@ module stdlib_linalg
552552
#:endfor
553553
end interface
554554

555+
interface eig
556+
!! Eigendecomposition of a square matrix: return eigenvalues, and optionally eigenvectors
557+
#:for rk,rt,ri in RC_KINDS_TYPES
558+
#:if rk!="xdp"
559+
module subroutine stdlib_linalg_eig_${ri}$(a,lambda,right,left,overwrite_a,err)
560+
!! Eigendecomposition of matrix A returning an array `lambda` of eigenvalues,
561+
!! and optionally right or left eigenvectors.
562+
!> Input matrix A[m,n]
563+
${rt}$, intent(inout), target :: a(:,:)
564+
!> Array of eigenvalues
565+
complex(${rk}$), intent(out) :: lambda(:)
566+
!> The columns of RIGHT contain the right eigenvectors of A
567+
complex(${rk}$), optional, intent(out), target :: right(:,:)
568+
!> The columns of LEFT contain the left eigenvectors of A
569+
complex(${rk}$), optional, intent(out), target :: left(:,:)
570+
!> [optional] Can A data be overwritten and destroyed?
571+
logical(lk), optional, intent(in) :: overwrite_a
572+
!> [optional] state return flag. On error if not requested, the code will stop
573+
type(linalg_state_type), optional, intent(out) :: err
574+
end subroutine stdlib_linalg_eig_${ri}$
575+
#:endif
576+
#:endfor
577+
end interface eig
578+
579+
interface eigvals
580+
!! Eigenvalues of a square matrix
581+
#:for rk,rt,ri in RC_KINDS_TYPES
582+
#:if rk!="xdp"
583+
module function stdlib_linalg_eigvals_${ri}$(a,err) result(lambda)
584+
!! Return an array of eigenvalues of matrix A.
585+
!> Input matrix A[m,n]
586+
${rt}$, intent(in), target :: a(:,:)
587+
!> [optional] state return flag. On error if not requested, the code will stop
588+
type(linalg_state_type), intent(out) :: err
589+
!> Array of singular values
590+
complex(${rk}$), allocatable :: lambda(:)
591+
end function stdlib_linalg_eigvals_${ri}$
592+
593+
module function stdlib_linalg_eigvals_noerr_${ri}$(a) result(lambda)
594+
!! Return an array of eigenvalues of matrix A.
595+
!> Input matrix A[m,n]
596+
${rt}$, intent(in), target :: a(:,:)
597+
!> Array of singular values
598+
complex(${rk}$), allocatable :: lambda(:)
599+
end function stdlib_linalg_eigvals_noerr_${ri}$
600+
#:endif
601+
#:endfor
602+
end interface eigvals
603+
604+
interface eigh
605+
!! Eigendecomposition of a real symmetric or complex hermitian matrix
606+
#:for rk,rt,ri in RC_KINDS_TYPES
607+
#:if rk!="xdp"
608+
module subroutine stdlib_linalg_eigh_${ri}$(a,lambda,vectors,upper_a,overwrite_a,err)
609+
!! Eigendecomposition of a real symmetric or complex Hermitian matrix A returning an array `lambda`
610+
!! of eigenvalues, and optionally right or left eigenvectors.
611+
!> Input matrix A[m,n]
612+
${rt}$, intent(inout), target :: a(:,:)
613+
!> Array of eigenvalues
614+
real(${rk}$), intent(out) :: lambda(:)
615+
!> The columns of vectors contain the orthonormal eigenvectors of A
616+
${rt}$, optional, intent(out), target :: vectors(:,:)
617+
!> [optional] Can A data be overwritten and destroyed?
618+
logical(lk), optional, intent(in) :: overwrite_a
619+
!> [optional] Should the upper/lower half of A be used? Default: lower
620+
logical(lk), optional, intent(in) :: upper_a
621+
!> [optional] state return flag. On error if not requested, the code will stop
622+
type(linalg_state_type), optional, intent(out) :: err
623+
end subroutine stdlib_linalg_eigh_${ri}$
624+
#:endif
625+
#:endfor
626+
end interface eigh
627+
628+
interface eigvalsh
629+
!! Eigenvalues of a real symmetric or complex hermitian matrix
630+
#:for rk,rt,ri in RC_KINDS_TYPES
631+
#:if rk!="xdp"
632+
module function stdlib_linalg_eigvalsh_${ri}$(a,upper_a,err) result(lambda)
633+
!! Return an array of eigenvalues of real symmetric / complex hermitian A
634+
!> Input matrix A[m,n]
635+
${rt}$, intent(in), target :: a(:,:)
636+
!> [optional] Should the upper/lower half of A be used? Default: lower
637+
logical(lk), optional, intent(in) :: upper_a
638+
!> [optional] state return flag. On error if not requested, the code will stop
639+
type(linalg_state_type), intent(out) :: err
640+
!> Array of singular values
641+
real(${rk}$), allocatable :: lambda(:)
642+
end function stdlib_linalg_eigvalsh_${ri}$
643+
644+
module function stdlib_linalg_eigvalsh_noerr_${ri}$(a,upper_a) result(lambda)
645+
!! Return an array of eigenvalues of real symmetric / complex hermitian A
646+
!> Input matrix A[m,n]
647+
${rt}$, intent(in), target :: a(:,:)
648+
!> [optional] Should the upper/lower half of A be used? Default: lower
649+
logical(lk), optional, intent(in) :: upper_a
650+
!> Array of singular values
651+
real(${rk}$), allocatable :: lambda(:)
652+
end function stdlib_linalg_eigvalsh_noerr_${ri}$
653+
#:endif
654+
#:endfor
655+
end interface eigvalsh
656+
555657
contains
556658

557659

src/stdlib_linalg_eigenvalues.fypp

Lines changed: 10 additions & 64 deletions
Original file line numberDiff line numberDiff line change
@@ -1,70 +1,15 @@
11
#:include "common.fypp"
22
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
3-
module stdlib_linalg_eigenvalues
3+
submodule (stdlib_linalg) stdlib_linalg_eigenvalues
44
!! Compute eigenvalues and eigenvectors
55
use stdlib_linalg_constants
6-
use stdlib_linalg_blas
7-
use stdlib_linalg_lapack
8-
use stdlib_linalg_state
9-
use iso_fortran_env,only:real32,real64,real128,int8,int16,int32,int64,stderr => error_unit
6+
use stdlib_linalg_lapack, only: geev, heev, syev
7+
use stdlib_linalg_state, only: linalg_state_type, linalg_error_handling, LINALG_ERROR, &
8+
LINALG_INTERNAL_ERROR, LINALG_VALUE_ERROR, LINALG_SUCCESS
109
implicit none(type,external)
11-
private
12-
13-
!> Eigendecomposition of a square matrix: return eigenvalues, and optionally eigenvectors
14-
public :: eig
15-
!> Eigenvalues of a square matrix
16-
public :: eigvals
17-
!> Eigendecomposition of a real symmetric or complex hermitian matrix
18-
public :: eigh
19-
!> Eigenvalues of a real symmetric or complex hermitian matrix
20-
public :: eigvalsh
21-
22-
! Numpy: eigenvalues, eigenvectors = eig(a)
23-
! eigenvalues = eigvals(a)
24-
! Scipy: eig(a, b=None, left=False, right=True, overwrite_a=False, overwrite_b=False, check_finite=True, homogeneous_eigvals=False)
25-
26-
! Numpy: eigenvalues, eigenvectors = eigh(a, uplo='L')
27-
! eigenvalues = eigvalsh(a)
28-
! Scipy: eigh(a, b=None, *, lower=True, eigvals_only=False, overwrite_a=False, overwrite_b=False, turbo=<object object>, eigvals=<object object>, type=1, check_finite=True, subset_by_index=None, subset_by_value=None, driver=None)
29-
30-
interface eig
31-
#:for rk,rt,ri in RC_KINDS_TYPES
32-
#:if rk!="xdp"
33-
module procedure stdlib_linalg_eig_${ri}$
34-
#:endif
35-
#:endfor
36-
end interface eig
37-
38-
interface eigvals
39-
#:for rk,rt,ri in RC_KINDS_TYPES
40-
#:if rk!="xdp"
41-
module procedure stdlib_linalg_eigvals_${ri}$
42-
module procedure stdlib_linalg_eigvals_noerr_${ri}$
43-
#:endif
44-
#:endfor
45-
end interface eigvals
4610

47-
interface eigh
48-
#:for rk,rt,ri in RC_KINDS_TYPES
49-
#:if rk!="xdp"
50-
module procedure stdlib_linalg_eigh_${ri}$
51-
#:endif
52-
#:endfor
53-
end interface eigh
54-
55-
interface eigvalsh
56-
#:for rk,rt,ri in RC_KINDS_TYPES
57-
#:if rk!="xdp"
58-
module procedure stdlib_linalg_eigvalsh_${ri}$
59-
module procedure stdlib_linalg_eigvalsh_noerr_${ri}$
60-
#:endif
61-
#:endfor
62-
end interface eigvalsh
63-
64-
6511
character(*), parameter :: this = 'eigenvalues'
6612

67-
6813
contains
6914

7015
!> Request for eigenvector calculation
@@ -236,13 +181,14 @@ module stdlib_linalg_eigenvalues
236181
lda = m
237182

238183
if (.not.(k>0 .and. m==n)) then
239-
err0 = linalg_state_type(this,LINALG_VALUE_ERROR,'invalid or matrix size a=',[m,n], &
240-
', must be square.')
184+
err0 = linalg_state_type(this,LINALG_VALUE_ERROR,&
185+
'invalid or matrix size a=',[m,n],', must be square.')
241186
call linalg_error_handling(err0,err)
242187
return
243188
elseif (.not.neig>=k) then
244-
err0 = linalg_state_type(this,LINALG_VALUE_ERROR,'eigenvalue array has insufficient size:',&
245-
' lambda=',neig,', n=',n)
189+
err0 = linalg_state_type(this,LINALG_VALUE_ERROR,&
190+
'eigenvalue array has insufficient size:',&
191+
' lambda=',neig,', n=',n)
246192
call linalg_error_handling(err0,err)
247193
return
248194
endif
@@ -590,4 +536,4 @@ module stdlib_linalg_eigenvalues
590536
#:endfor
591537

592538

593-
end module stdlib_linalg_eigenvalues
539+
end submodule stdlib_linalg_eigenvalues

0 commit comments

Comments
 (0)