The Elk Code
 
Loading...
Searching...
No Matches
eveqnfvz.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 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 eveqnfvz(nmatp,h,o,evalfv,evecfv)
7use modmain
8use modomp
9implicit none
10! arguments
11integer, intent(in) :: nmatp
12complex(8), intent(in) :: h(*),o(*)
13real(8), intent(out) :: evalfv(nstfv)
14complex(8), intent(out) :: evecfv(nmatmax,nstfv)
15! local variables
16integer nb,lwork,nts
17integer m,info,nthd
18real(8) vl,vu
19real(8) ts0,ts1
20! automatic arrays
21integer iwork(5*nmatp),ifail(nmatp)
22real(8) w(nmatp),rwork(7*nmatp)
23! allocatable arrays
24complex(8), allocatable :: work(:)
25! external functions
26integer, external :: ilaenv
27call timesec(ts0)
28! find the optimal blocksize for allocating the work array
29nb=ilaenv(1,'ZHETRD','U',nmatp,-1,-1,-1)
30nb=max(nb,1)
31lwork=(nb+1)*nmatp
32allocate(work(lwork))
33! enable MKL parallelism
34call holdthd(maxthdmkl,nthd)
35nts=mkl_set_num_threads_local(nthd)
36! diagonalise the matrix
37call zhegvx(1,'V','I','U',nmatp,h,nmatp,o,nmatp,vl,vu,1,nstfv,evaltol,m,w, &
38 evecfv,nmatmax,work,lwork,rwork,iwork,ifail,info)
39nts=mkl_set_num_threads_local(0)
40call freethd(nthd)
41if (info /= 0) then
42 write(*,*)
43 write(*,'("Error(eveqnfvz): diagonalisation failed")')
44 write(*,'(" ZHEGVX returned INFO = ",I8)') info
45 if (info > nmatp) then
46 write(*,'(" The leading minor of the overlap matrix of order ",I8)') &
47 info-nmatp
48 write(*,'(" is not positive definite")')
49 write(*,'(" Order of overlap matrix : ",I8)') nmatp
50 end if
51 write(*,*)
52 stop
53end if
54evalfv(1:nstfv)=w(1:nstfv)
55deallocate(work)
56call timesec(ts1)
57!$OMP ATOMIC
58timefv=timefv+ts1-ts0
59end subroutine
60
subroutine eveqnfvz(nmatp, h, o, evalfv, evecfv)
Definition eveqnfvz.f90:7
real(8) evaltol
Definition modmain.f90:916
real(8) timefv
Definition modmain.f90:1216
integer maxthdmkl
Definition modomp.f90:15
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine timesec(ts)
Definition timesec.f90:10