@@ -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+
555657contains
556658
557659
0 commit comments