The Elk Code
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 
6 subroutine eveqnfvz(nmatp,h,o,evalfv,evecfv)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: nmatp
12 complex(8), intent(in) :: h(*),o(*)
13 real(8), intent(out) :: evalfv(nstfv)
14 complex(8), intent(out) :: evecfv(nmatmax,nstfv)
15 ! local variables
16 integer nb,lwork,nts
17 integer m,info,nthd
18 real(8) vl,vu
19 real(8) ts0,ts1
20 ! automatic arrays
21 integer iwork(5*nmatp),ifail(nmatp)
22 real(8) w(nmatp),rwork(7*nmatp)
23 ! allocatable arrays
24 complex(8), allocatable :: work(:)
25 ! external functions
26 integer, external :: ilaenv
27 call timesec(ts0)
28 ! find the optimal blocksize for allocating the work array
29 nb=ilaenv(1,'ZHETRD','U',nmatp,-1,-1,-1)
30 nb=max(nb,1)
31 lwork=(nb+1)*nmatp
32 allocate(work(lwork))
33 ! enable MKL parallelism
34 call holdthd(maxthdmkl,nthd)
35 nts=mkl_set_num_threads_local(nthd)
36 ! diagonalise the matrix
37 call 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)
39 nts=mkl_set_num_threads_local(0)
40 call freethd(nthd)
41 if (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
53 end if
54 evalfv(1:nstfv)=w(1:nstfv)
55 deallocate(work)
56 call timesec(ts1)
57 !$OMP ATOMIC
58 timefv=timefv+ts1-ts0
59 end subroutine
60 
subroutine eveqnfvz(nmatp, h, o, evalfv, evecfv)
Definition: eveqnfvz.f90:7
Definition: modomp.f90:6
integer maxthdmkl
Definition: modomp.f90:15
real(8) evaltol
Definition: modmain.f90:919
subroutine timesec(ts)
Definition: timesec.f90:10
real(8) timefv
Definition: modmain.f90:1219
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78