The Elk Code
 
Loading...
Searching...
No Matches
exxengyk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine exxengyk(ikp)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ikp
11! local variables
12integer iq,ik,jk,m
13integer nst1,nst2,ist,jst
14integer is,ia,ias
15integer nrc,nrci,npc
16integer iv(3),ig
17real(8) ex,vc(3)
18complex(8) z1
19! automatic arrays
20integer idx(nstsv)
21! allocatable arrays
22real(8), allocatable :: vgqc(:,:),gqc(:),gclgq(:),jlgqrmt(:,:,:)
23complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
24complex(8), allocatable :: ylmgq(:,:),sfacgq(:,:)
25complex(4), allocatable :: wfmt1(:,:,:,:),wfir1(:,:,:),wfcr(:,:)
26complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
27complex(4), allocatable :: crhomt(:,:),crhoir(:)
28complex(4), allocatable :: cvclmt(:,:),cvclir(:)
29! external functions
30complex(8), external :: zcfinp,zcfmtinp
31! get the eigenvectors from file for input reduced k-point
32allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
33call getevecfv(filext,ikp,vkl(:,ikp),vgkl(:,:,:,ikp),evecfv)
34call getevecsv(filext,ikp,vkl(:,ikp),evecsv)
35! find the matching coefficients
36allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
37call match(ngk(1,ikp),vgkc(:,:,1,ikp),gkc(:,1,ikp),sfacgk(:,:,1,ikp),apwalm)
38! count and index the occupied states
39nst1=0
40do ist=1,nstsv
41 if (evalsv(ist,ikp) > efermi) cycle
42 nst1=nst1+1
43 idx(nst1)=ist
44end do
45! calculate the wavefunctions for occupied states of the input k-point
46allocate(wfmt1(npcmtmax,natmtot,nspinor,nst1),wfir1(ngtc,nspinor,nst1))
47call genwfsv_sp(.false.,.false.,nst1,idx,ngdgc,igfc,ngk(1,ikp),igkig(:,1,ikp), &
48 apwalm,evecfv,evecsv,wfmt1,ngtc,wfir1)
49! allocate local arrays
50allocate(vgqc(3,ngvc),gqc(ngvc),gclgq(ngvc))
51allocate(jlgqrmt(0:lnpsd,ngvc,nspecies))
52allocate(ylmgq(lmmaxo,ngvc),sfacgq(ngvc,natmtot))
53allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv))
54allocate(wfir2(ngtc,nspinor,nstsv))
55allocate(crhomt(npcmtmax,natmtot),crhoir(ngtc))
56allocate(cvclmt(npcmtmax,natmtot),cvclir(ngtc))
57! zero the local exchange energy variable
58ex=0.d0
59! start loop over non-reduced k-point set
60do ik=1,nkptnr
61! equivalent reduced k-point
62 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
63! determine the q-vector
64 iv(:)=ivk(:,ikp)-ivk(:,ik)
65 iv(:)=modulo(iv(:),ngridk(:))
66! check if the q-point is in user-defined set
67 iv(:)=iv(:)*ngridq(:)
68 if (any(mod(iv(:),ngridk(:)) /= 0)) cycle
69 iv(:)=iv(:)/ngridk(:)
70 iq=ivqiq(iv(1),iv(2),iv(3))
71 vc(:)=vkc(:,ikp)-vkc(:,ik)
72 do ig=1,ngvc
73! determine the G+q-vectors
74 vgqc(1:3,ig)=vgc(1:3,ig)+vc(1:3)
75! G+q-vector length
76 gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
77! spherical harmonics for G+q-vectors
78 call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
79 end do
80! structure factors for G+q
81 call gensfacgp(ngvc,vgqc,ngvc,sfacgq)
82! generate the regularised Coulomb Green's function in G+q-space
83 call gengclgq(.true.,iq,ngvc,gqc,gclgq)
84! compute the required spherical Bessel functions
85 call genjlgprmt(lnpsd,ngvc,gqc,ngvc,jlgqrmt)
86! find the matching coefficients
87 call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
88! get the eigenvectors from file for non-reduced k-point
89 call getevecfv(filext,0,vkl(:,ik),vgkl(:,:,1,ik),evecfv)
90 call getevecsv(filext,0,vkl(:,ik),evecsv)
91! count and index the occupied states
92 nst2=0
93 do jst=1,nstsv
94 if (evalsv(jst,jk) > efermi) cycle
95 nst2=nst2+1
96 idx(nst2)=jst
97 end do
98! calculate the wavefunctions for occupied states
99 call genwfsv_sp(.false.,.false.,nst2,idx,ngdgc,igfc,ngk(1,ik),igkig(:,1,ik), &
100 apwalm,evecfv,evecsv,wfmt2,ngtc,wfir2)
101!--------------------------------------------!
102! valence-valence-valence contribution !
103!--------------------------------------------!
104 do jst=1,nst2
105 do ist=1,nst1
106! calculate the complex overlap density
107 call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,jst),wfir2(:,:,jst), &
108 wfmt1(:,:,:,ist),wfir1(:,:,ist),crhomt,crhoir)
109! calculate the Coulomb potential
112 gclgq,ngvc,jlgqrmt,ylmgq,sfacgq,crhoir,npcmtmax,cvclmt,cvclir)
113 cvclir(1:ngtc)=cvclir(1:ngtc)*cfrc(1:ngtc)
114 z1=zcfinp(crhomt,crhoir,cvclmt,cvclir)
115 ex=ex-0.5d0*occmax*wkpt(ikp)*wqptnr*z1%re
116 end do
117 end do
118! end loop over non-reduced k-point set
119end do
120deallocate(vgqc,gqc,gclgq,jlgqrmt)
121deallocate(evecfv,evecsv)
122deallocate(apwalm,ylmgq,sfacgq)
123deallocate(wfmt2,wfir2)
124!-----------------------------------------!
125! valence-core-valence contribution !
126!-----------------------------------------!
127allocate(wfcr(npcmtmax,2))
128! begin loops over atoms and species
129do is=1,nspecies
130 nrc=nrcmt(is)
131 nrci=nrcmti(is)
132 npc=npcmt(is)
133 do ia=1,natoms(is)
134 ias=idxas(ia,is)
135 do jst=1,nstsp(is)
136 if (spcore(jst,is)) then
137 do m=-ksp(jst,is),ksp(jst,is)-1
138! generate the core wavefunction in spherical coordinates (pass in m-1/2)
139 call wavefcr(.false.,lradstp,is,ia,jst,m,npcmtmax,wfcr)
140 do ist=1,nst1
141! calculate the complex overlap density in spherical harmonics
142 if (spinpol) then
143 call crho2(npc,wfcr,wfcr(:,2),wfmt1(:,ias,1,ist), &
144 wfmt1(:,ias,2,ist),crhomt(:,ias))
145 else
146 call crho1(npc,wfcr,wfmt1(:,ias,1,ist),crhomt(:,ias))
147 end if
148 call cfshtip(nrc,nrci,crhomt(:,ias))
149! calculate the Coulomb potential
150 call cpotclmt(nrc,nrci,nrcmtmax,rlcmt(:,:,is),wprcmt(:,:,is), &
151 crhomt(:,ias),cvclmt(:,ias))
152 z1=zcfmtinp(nrc,nrci,wr2cmt(:,is),crhomt(:,ias),cvclmt(:,ias))
153 ex=ex-occmax*wkpt(ikp)*z1%re
154 end do
155! end loop over m
156 end do
157! end loop over jst
158 end if
159 end do
160! end loops over atoms and species
161 end do
162end do
163! add to global exchange energy
164!$OMP CRITICAL(exxengyk_)
165engyx=engyx+ex
166!$OMP END CRITICAL(exxengyk_)
167deallocate(wfmt1,wfir1,wfcr)
168deallocate(crhomt,crhoir,cvclmt,cvclir)
169return
170
171contains
172
173pure subroutine crho1(n,wf1,wf2,crho)
174implicit none
175integer, intent(in) :: n
176complex(4), intent(in) :: wf1(n),wf2(n)
177complex(4), intent(out) :: crho(n)
178crho(1:n)=conjg(wf1(1:n))*wf2(1:n)
179end subroutine
180
181pure subroutine crho2(n,wf11,wf12,wf21,wf22,crho)
182implicit none
183integer, intent(in) :: n
184complex(4), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
185complex(4), intent(out) :: crho(n)
186crho(1:n)=conjg(wf11(1:n))*wf21(1:n)+conjg(wf12(1:n))*wf22(1:n)
187end subroutine
188
189end subroutine
190
subroutine cfshtip(nr, nri, cfmt)
Definition cfshtip.f90:7
pure subroutine cpotclmt(nr, nri, ld, rl, wpr, crhomt, cvclmt)
Definition cpotclmt.f90:7
subroutine cpotcoul(nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, crhoir, ld3, cvclmt, cvclir)
Definition cpotcoul.f90:8
pure subroutine crho2(n, wf11, wf12, wf21, wf22, crho)
Definition exxengy.f90:86
pure subroutine crho1(n, wf1, wf2, crho)
Definition exxengyk.f90:174
subroutine exxengyk(ikp)
Definition exxengyk.f90:7
subroutine gencrho(tsh, tspc, ngt, wfmt1, wfir1, wfmt2, wfir2, crhomt, crhoir)
Definition gencrho.f90:7
subroutine gencvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, crhomt, cvclmt)
Definition gencvclmt.f90:7
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition gengclgq.f90:7
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition gensfacgp.f90:10
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv_sp.f90:8
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
real(8) efermi
Definition modmain.f90:904
integer nkptnr
Definition modmain.f90:463
integer nspinor
Definition modmain.f90:267
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
character(256) filext
Definition modmain.f90:1300
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
logical spinpol
Definition modmain.f90:228
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer lradstp
Definition modmain.f90:171
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
integer lmmaxapw
Definition modmain.f90:199
integer apwordmax
Definition modmain.f90:760
integer ngtc
Definition modmain.f90:392
real(8), dimension(:,:,:), allocatable rlcmt
Definition modmain.f90:181
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer npcmtmax
Definition modmain.f90:216
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(:,:,:), allocatable wprcmt
Definition modmain.f90:191
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer nmatmax
Definition modmain.f90:848
real(8) engyx
Definition modmain.f90:972
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
integer lmaxo
Definition modmain.f90:201
real(8) wqptnr
Definition modmain.f90:551
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
real(8) occmax
Definition modmain.f90:898
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer nrcmtmax
Definition modmain.f90:175
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
integer lmmaxo
Definition modmain.f90:203
integer ngvc
Definition modmain.f90:398
integer nstfv
Definition modmain.f90:884
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
integer, dimension(3) ngridq
Definition modmain.f90:515
integer lnpsd
Definition modmain.f90:628
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
pure subroutine wavefcr(tsh, lrstp, is, ia, ist, m, ld, wfcr)
Definition wavefcr.f90:7
pure complex(8) function zcfmtinp(nr, nri, wr, cfmt1, cfmt2)
Definition zcfmtinp.f90:7