The Elk Code
eveqnhf.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2006 J. K. Dewhurst and S. Sharma.
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 eveqnhf(ikp,vmt,vir,bmt,bir,evecsvp)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: ikp
12 real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtc)
13 real(8), intent(in) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
14 complex(8), intent(inout) :: evecsvp(nstsv,nstsv)
15 ! local variables
16 integer ik,jk,nst
17 integer ist1,ist2,ist3,jst3
18 integer iv(3),iq,ig,nthd
19 real(8) vc(3),t1
20 complex(8) z1
21 ! automatic arrays
22 integer idx(nstsv)
23 real(8) vgqc(3,ngvc),gqc(ngvc),gclgq(ngvc)
24 ! allocatable arrays
25 real(8), allocatable :: jlgqrmt(:,:,:)
26 complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
27 complex(8), allocatable :: ylmgq(:,:),sfacgq(:,:)
28 complex(8), allocatable :: h(:,:),v(:,:),kmat(:,:)
29 complex(4), allocatable :: wfmt1(:,:,:,:),wfir1(:,:,:)
30 complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
31 complex(4), allocatable :: crhomt(:,:,:),crhoir(:,:)
32 complex(4), allocatable :: cvclmt(:,:),cvclir(:)
33 ! external functions
34 complex(8), external :: zcfinp
35 !$OMP CRITICAL(eveqnhf_)
36 write(*,'("Info(eveqnhf): ",I0," of ",I0," k-points")') ikp,nkpt
37 !$OMP END CRITICAL(eveqnhf_)
38 ! allocate local arrays
39 allocate(jlgqrmt(0:lnpsd,ngvc,nspecies))
40 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
41 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
42 allocate(ylmgq(lmmaxo,ngvc),sfacgq(ngvc,natmtot))
43 allocate(h(nstsv,nstsv),v(nstsv,nstsv))
44 allocate(wfmt1(npcmtmax,natmtot,nspinor,nstsv),wfir1(ngtc,nspinor,nstsv))
45 allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv),wfir2(ngtc,nspinor,nstsv))
46 allocate(crhomt(npcmtmax,natmtot,nstsv),crhoir(ngtc,nstsv))
47 ! get the first-variational eigenvectors from file for input reduced k-point
48 call getevecfv(filext,ikp,vkl(:,ikp),vgkl(:,:,1,ikp),evecfv)
49 ! find the matching coefficients
50 call match(ngk(1,ikp),vgkc(:,:,1,ikp),gkc(:,1,ikp),sfacgk(:,:,1,ikp),apwalm)
51 ! calculate the wavefunctions for all states of the input k-point
52 call genwfsv_sp(.false.,.true.,nstsv,[0],ngdgc,igfc,ngk(1,ikp),igkig(:,1,ikp), &
53  apwalm,evecfv,evecsvp,wfmt1,ngtc,wfir1)
54 !-----------------------------------------!
55 ! local potential matrix elements !
56 !-----------------------------------------!
57 if (hybrid.and.spinpol) then
58 ! magnetic field matrix elements in hybrid case
59  call genvbmatk(vmt,vir,bmt,bir,ngk(1,ikp),igkig(:,1,ikp),wfmt1,ngtc,wfir1,h)
60 else
61  call genvmatk(vmt,vir,ngk(1,ikp),igkig(:,1,ikp),wfmt1,ngtc,wfir1,h)
62 end if
63 ! Fourier transform wavefunctions to real-space
64 call cftwfir(ngk(1,ikp),igkig(:,1,ikp),wfir1)
65 !---------------------------------!
66 ! kinetic matrix elements !
67 !---------------------------------!
68 allocate(kmat(nstsv,nstsv))
69 call getkmat(ikp,kmat)
70 call zgemm('N','N',nstsv,nstsv,nstsv,zone,kmat,nstsv,evecsvp,nstsv,zzero,v, &
71  nstsv)
72 call zgemm('C','N',nstsv,nstsv,nstsv,zone,evecsvp,nstsv,v,nstsv,zone,h,nstsv)
73 deallocate(kmat)
74 !------------------------------!
75 ! Fock matrix elements !
76 !------------------------------!
77 v(:,:)=0.d0
78 ! loop over non-reduced k-point set
79 do ik=1,nkptnr
80 ! equivalent reduced k-point
81  jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
82 ! determine the q-vector
83  iv(1:3)=ivk(1:3,ikp)-ivk(1:3,ik)
84  iv(1:3)=modulo(iv(1:3),ngridk(1:3))
85 ! check if the q-point is in user-defined set
86  iv(1:3)=iv(1:3)*ngridq(1:3)
87  if (any(mod(iv(1:3),ngridk(1:3)) /= 0)) cycle
88  iv(1:3)=iv(1:3)/ngridk(1:3)
89  iq=ivqiq(iv(1),iv(2),iv(3))
90  vc(1:3)=vkc(1:3,ikp)-vkc(1:3,ik)
91  do ig=1,ngvc
92 ! determine the G+q-vectors
93  vgqc(1:3,ig)=vgc(1:3,ig)+vc(1:3)
94 ! G+q-vector length
95  gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
96 ! spherical harmonics for G+q-vectors
97  call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
98  end do
99 ! structure factors for G+q
100  call gensfacgp(ngvc,vgqc,ngvc,sfacgq)
101 ! generate the regularised Coulomb Green's function in G+q-space
102  call gengclgq(.true.,iq,ngvc,gqc,gclgq)
103 ! compute the required spherical Bessel functions
104  call genjlgprmt(lnpsd,ngvc,gqc,ngvc,jlgqrmt)
105 ! find the matching coefficients
106  call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
107 ! get the eigenvectors from file for non-reduced k-point
108  call getevecfv(filext,0,vkl(:,ik),vgkl(:,:,1,ik),evecfv)
109  call getevecsv(filext,0,vkl(:,ik),evecsv)
110 ! count and index occupied states
111  nst=0
112  do ist3=1,nstsv
113  if (abs(occsv(ist3,jk)) < epsocc) cycle
114  nst=nst+1
115  idx(nst)=ist3
116  end do
117 ! calculate the wavefunctions for occupied states
118  call genwfsv_sp(.false.,.false.,nst,idx,ngdgc,igfc,ngk(1,ik),igkig(:,1,ik), &
119  apwalm,evecfv,evecsv,wfmt2,ngtc,wfir2)
120  call holdthd(nstsv,nthd)
121 !$OMP PARALLEL DEFAULT(SHARED) &
122 !$OMP PRIVATE(cvclmt,cvclir) &
123 !$OMP PRIVATE(ist1,ist2,ist3,jst3,t1,z1) &
124 !$OMP NUM_THREADS(nthd)
125  allocate(cvclmt(npcmtmax,natmtot),cvclir(ngtc))
126  do ist3=1,nst
127  jst3=idx(ist3)
128 ! calculate the complex overlap densities for all states (T. McQueen)
129 !$OMP DO SCHEDULE(DYNAMIC)
130  do ist1=1,nstsv
131  call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,ist3),wfir2(:,:,ist3), &
132  wfmt1(:,:,:,ist1),wfir1(:,:,ist1),crhomt(:,:,ist1),crhoir(:,ist1))
133  end do
134 !$OMP END DO
135  t1=wqptnr*occsv(jst3,jk)/occmax
136 !$OMP DO SCHEDULE(DYNAMIC)
137  do ist2=1,nstsv
138 ! calculate the Coulomb potential
139  call gencvclmt(nrcmt,nrcmti,nrcmtmax,rlcmt,wprcmt,npcmtmax, &
140  crhomt(:,:,ist2),cvclmt)
141  call cpotcoul(nrcmt,nrcmti,npcmt,nrcmtmax,rlcmt,ngdgc,igfc,ngvc,gqc, &
142  gclgq,ngvc,jlgqrmt,ylmgq,sfacgq,crhoir(:,ist2),npcmtmax,cvclmt,cvclir)
143  cvclir(:)=cvclir(:)*cfrc(:)
144  do ist1=1,ist2
145  z1=zcfinp(crhomt(:,:,ist1),crhoir(:,ist1),cvclmt,cvclir)
146  v(ist1,ist2)=v(ist1,ist2)-t1*z1
147  end do
148  end do
149 !$OMP END DO
150  end do
151  deallocate(cvclmt,cvclir)
152 !$OMP END PARALLEL
153  call freethd(nthd)
154 ! end loop over non-reduced k-point set
155 end do
156 deallocate(jlgqrmt,ylmgq,sfacgq,apwalm,evecfv)
157 deallocate(wfmt1,wfir1,wfmt2,wfir2,crhomt,crhoir)
158 ! scale the Coulomb matrix elements in the case of a hybrid functional
159 if (hybrid) v(:,:)=hybridc*v(:,:)
160 ! add the Coulomb matrix elements to Hamiltonian
161 h(:,:)=h(:,:)+v(:,:)
162 !----------------------------------------------!
163 ! diagonalise Hartree-Fock Hamiltonian !
164 !----------------------------------------------!
165 call eveqnzh(nstsv,nstsv,h,evalsv(:,ikp))
166 ! apply unitary transformation to the third-variational states so that they
167 ! refer to the first-variational basis
168 evecsv(:,:)=evecsvp(:,:)
169 call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,h,nstsv,zzero,evecsvp, &
170  nstsv)
171 deallocate(evecsv,h,v)
172 end subroutine
173 
integer nmatmax
Definition: modmain.f90:853
subroutine gencvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, crhomt, cvclmt)
Definition: gencvclmt.f90:7
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1299
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:919
integer lmmaxo
Definition: modmain.f90:203
logical spinpol
Definition: modmain.f90:228
integer lmmaxapw
Definition: modmain.f90:199
logical hybrid
Definition: modmain.f90:1150
integer nkpt
Definition: modmain.f90:461
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
real(8) wqptnr
Definition: modmain.f90:551
subroutine gencrho(tsh, tspc, ngt, wfmt1, wfir1, wfmt2, wfir2, crhomt, crhoir)
Definition: gencrho.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
Definition: modomp.f90:6
real(8) epsocc
Definition: modmain.f90:901
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
complex(8), parameter zone
Definition: modmain.f90:1238
subroutine getkmat(ik, kmat)
Definition: getkmat.f90:7
integer lmaxo
Definition: modmain.f90:201
integer nkptnr
Definition: modmain.f90:463
subroutine eveqnzh(n, ld, a, w)
Definition: eveqnzh.f90:7
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
real(8), dimension(:,:), allocatable vkc
Definition: modmain.f90:473
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer nrcmtmax
Definition: modmain.f90:175
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition: gengclgq.f90:7
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine genvmatk(vmt, vir, ngp, igpig, wfmt, ld, wfgp, vmat)
Definition: genvmatk.f90:7
real(8) occmax
Definition: modmain.f90:899
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv_sp.f90:8
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
real(8), dimension(:,:,:), allocatable rlcmt
Definition: modmain.f90:181
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:903
real(8), dimension(:,:,:), allocatable wprcmt
Definition: modmain.f90:191
integer nspinor
Definition: modmain.f90:267
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
complex(8), parameter zzero
Definition: modmain.f90:1238
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine cpotcoul(nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, crhoir, ld3, cvclmt, cvclir)
Definition: cpotcoul.f90:8
integer apwordmax
Definition: modmain.f90:760
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
integer lnpsd
Definition: modmain.f90:628
integer nspecies
Definition: modmain.f90:34
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
Definition: genjlgprmt.f90:10
subroutine cftwfir(ngp, igpig, wfir)
Definition: cftwfir.f90:7
real(8), dimension(:), allocatable cfrc
Definition: modmain.f90:438
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
subroutine genvbmatk(vmt, vir, bmt, bir, ngp, igpig, wfmt, ld, wfgp, vbmat)
Definition: genvbmatk.f90:7
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:885
real(8) hybridc
Definition: modmain.f90:1152
subroutine eveqnhf(ikp, vmt, vir, bmt, bir, evecsvp)
Definition: eveqnhf.f90:7
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465