@brief Find eigenvalues and eigenvectors of symmetric matrix in packed format @param[in] mode algorithm of diagonalization (not used now) @param[in] n matrix dimension @param[in] ldvect leading dimension of eigenvector matrix @param[in] nvect required number of eigenvectors @param[in,out] h matrix to be diagonalized @param[out] eig eigenvalues @param[out] vector eigenvectors @param[out] ierr status
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | mode | |||
integer, | intent(in) | :: | ldvect | |||
integer, | intent(in) | :: | nvect | |||
integer, | intent(in) | :: | n | |||
real(kind=dp), | intent(inout) | :: | h(*) | |||
real(kind=dp), | intent(out) | :: | eig(*) | |||
real(kind=dp), | intent(out) | :: | vector(*) | |||
integer, | intent(out), | optional | :: | ierr |
subroutine diag_symm_packed(mode, ldvect, nvect, n, h, eig, vector, ierr) use messages, only: show_message, WITH_ABORT, WITHOUT_ABORT ! integer, intent(in) :: mode integer, intent(in) :: ldvect, nvect, n integer, optional, intent(out) :: ierr real(dp), intent(inout) :: h(*) real(kind=dp), intent(out) :: eig(*), vector(*) integer(blas_int), dimension(:), allocatable :: iwork, ifail integer(blas_int) :: ldvect_, n_, nvect_, info, ione integer :: iok real(dp), dimension(:), allocatable :: work real(dp) :: abstol, dlamch logical :: fatal character(16) :: driver ldvect_ = int(ldvect, kind=blas_int) n_ = int(n, kind=blas_int) nvect_ = int(nvect, kind=blas_int) ione = 1 fatal = WITH_ABORT if (present(ierr)) fatal = WITHOUT_ABORT allocate (work(n*8), iwork(n*5), ifail(n), stat=iok) if (iok /= 0) then if (present(ierr)) ierr = iok call show_message('Cannot allocate memory', fatal) return end if if (nvect == n .and. ldvect >= n) then driver = 'DSPEV' call dspev('V', 'U', n_, h, eig, vector, ldvect_, work, info) else abstol = two*DLAMCH('S') driver = 'DSPEVX' call dspevx('V', 'A', 'U', & ldvect_, h, zero, zero, ione, ione, abstol, n_, & eig, vector, nvect_, work, iwork, ifail, info) end if if (present(ierr)) ierr = info if (info /= 0) then call show_message('(A,I0)', & trim(driver)//' FAILED! INFO: ', int(info), fatal) end if end subroutine diag_symm_packed