The Elk Code
eveqnzhg.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 eveqnzhg(n,m,ld1,a,b,w,ld2,z)
7 use modomp
8 implicit none
9 ! arguments
10 integer, intent(in) :: n,m,ld1
11 complex(8), intent(in) :: a(ld1,*),b(ld1,*)
12 real(8), intent(out) :: w(m)
13 integer, intent(in) :: ld2
14 complex(8), intent(out) :: z(ld2,m)
15 ! local variables
16 integer nb,lwork,nts
17 integer p,info,nthd
18 real(8) vl,vu
19 ! automatic arrays
20 integer iwork(5*n),ifail(n)
21 real(8) wn(n),rwork(7*n)
22 ! external functions
23 integer, external :: ilaenv
24 ! find the optimal blocksize for allocating the work array
25 nb=ilaenv(1,'ZHETRD','U',n,-1,-1,-1)
26 nb=max(nb,1)
27 lwork=(nb+1)*n
28 ! enable MKL parallelism
29 call holdthd(maxthdmkl,nthd)
30 nts=mkl_set_num_threads_local(nthd)
31 block
32 complex(8) work(lwork)
33 ! diagonalise the matrix
34 call zhegvx(1,'V','I','U',n,a,ld1,b,ld1,vl,vu,1,m,-1.d0,p,wn,z,ld2,work,lwork, &
35  rwork,iwork,ifail,info)
36 end block
37 nts=mkl_set_num_threads_local(0)
38 call freethd(nthd)
39 if (info /= 0) then
40  write(*,*)
41  write(*,'("Error(eveqnzhg): diagonalisation failed")')
42  write(*,'(" ZHEGVX returned INFO = ",I0)') info
43  if (info > n) then
44  write(*,'(" The leading minor of the overlap matrix of order ",I0)') info-n
45  write(*,'(" is not positive definite")')
46  write(*,'(" Order of overlap matrix : ",I0)') n
47  end if
48  write(*,*)
49  stop
50 end if
51 w(1:m)=wn(1:m)
52 end subroutine
53 
Definition: modomp.f90:6
integer maxthdmkl
Definition: modomp.f90:15
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine eveqnzhg(n, m, ld1, a, b, w, ld2, z)
Definition: eveqnzhg.f90:7