The Elk Code
 
Loading...
Searching...
No Matches
eveqnzh.f90
Go to the documentation of this file.
1
2! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine eveqnzh(n,ld,a,w)
7use modmain
8use modomp
9implicit none
10! arguments
11integer, intent(in) :: n,ld
12complex(8), intent(inout) :: a(ld,n)
13real(8), intent(out) :: w(n)
14! local variables
15integer nb,lwork,lrwork
16integer nts,info,nthd
17! automatic arrays
18integer iwork(3+5*n)
19! allocatable arrays
20real(8), allocatable :: rwork(:)
21complex(8), allocatable :: work(:)
22! external functions
23integer, external :: ilaenv
24! find the optimal blocksize for allocating the work array
25nb=ilaenv(1,'ZHETRD','U',n,-1,-1,-1)
26lwork=max(2*n+n**2,(nb+1)*n)
27lrwork=1+5*n+2*n**2
28allocate(work(lwork),rwork(lrwork))
29! enable MKL parallelism
30call holdthd(maxthdmkl,nthd)
31nts=mkl_set_num_threads_local(nthd)
32call zheevd('V','U',n,a,ld,w,work,lwork,rwork,lrwork,iwork,3+5*n,info)
33nts=mkl_set_num_threads_local(0)
34call freethd(nthd)
35if (info /= 0) then
36 write(*,*)
37 write(*,'("Error(eveqnzh): diagonalisation failed")')
38 write(*,'(" ZHEEVD returned INFO = ",I8)') info
39 write(*,*)
40 stop
41end if
42deallocate(rwork,work)
43end subroutine
44
subroutine eveqnzh(n, ld, a, w)
Definition eveqnzh.f90:7
integer maxthdmkl
Definition modomp.f90:15
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106