The Elk Code
 
Loading...
Searching...
No Matches
oepvclk.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 oepvclk(ikp,vclcv,vclvv)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ikp
11complex(8), intent(out) :: vclcv(ncrmax,natmtot,nstsv)
12complex(8), intent(out) :: vclvv(nstsv,nstsv)
13! local variables
14integer ik,jk,nst,ist1,ist2,ist3
15integer is,ia,ias,nrc,nrci,npc
16integer iv(3),ig,iq,ic,jc,m1,m2
17real(8) 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(:,:,:)
26complex(4), allocatable :: wfmt2(:,:,:,:),wfir2(:,:,:)
27complex(4), allocatable :: wfcr1(:,:),wfcr2(:,:)
28complex(4), allocatable :: crhomt1(:,:,:),crhomt2(:,:),crhoir1(:,:)
29complex(4), allocatable :: cvclmt(:,:),cvclir(:)
30! external functions
31complex(8), external :: zcfinp,zcfmtinp
32! allocate local arrays
33allocate(vgqc(3,ngvc),gqc(ngvc),gclgq(ngvc))
34allocate(jlgqrmt(0:lnpsd,ngvc,nspecies))
35allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
36allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
37allocate(ylmgq(lmmaxo,ngvc),sfacgq(ngvc,natmtot))
38allocate(wfmt1(npcmtmax,natmtot,nspinor,nstsv),wfir1(ngtc,nspinor,nstsv))
39allocate(wfmt2(npcmtmax,natmtot,nspinor,nstsv),wfir2(ngtc,nspinor,nstsv))
40allocate(wfcr1(npcmtmax,2),wfcr2(npcmtmax,2))
41allocate(crhomt1(npcmtmax,natmtot,nstsv),crhoir1(ngtc,nstsv))
42allocate(crhomt2(npcmtmax,nstcr))
43allocate(cvclmt(npcmtmax,natmtot),cvclir(ngtc))
44! zero the Coulomb matrix elements
45vclcv(:,:,:)=0.d0
46vclvv(:,:)=0.d0
47! get the eigenvectors from file for input reduced k-point
48call getevecfv(filext,ikp,vkl(:,ikp),vgkl(:,:,:,ikp),evecfv)
49call getevecsv(filext,ikp,vkl(:,ikp),evecsv)
50! find the matching coefficients
51call match(ngk(1,ikp),vgkc(:,:,1,ikp),gkc(:,1,ikp),sfacgk(:,:,1,ikp),apwalm)
52! calculate the wavefunctions for all states of the input k-point
53call genwfsv_sp(.false.,.false.,nstsv,[0],ngdgc,igfc,ngk(1,ikp),igkig(:,1,ikp),&
54 apwalm,evecfv,evecsv,wfmt1,ngtc,wfir1)
55! loop over non-reduced k-point set
56do ik=1,nkptnr
57! equivalent reduced k-point
58 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
59! determine the q-vector
60 iv(:)=ivk(:,ikp)-ivk(:,ik)
61 iv(:)=modulo(iv(:),ngridk(:))
62! check if the q-point is in user-defined set
63 iv(:)=iv(:)*ngridq(:)
64 if (any(mod(iv(:),ngridk(:)) /= 0)) cycle
65 iv(:)=iv(:)/ngridk(:)
66 iq=ivqiq(iv(1),iv(2),iv(3))
67 vc(:)=vkc(:,ikp)-vkc(:,ik)
68 do ig=1,ngvc
69! determine the G+q-vectors
70 vgqc(1:3,ig)=vgc(1:3,ig)+vc(1:3)
71! G+q-vector length
72 gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
73! spherical harmonics for G+q-vectors
74 call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
75 end do
76! structure factors for G+q
77 call gensfacgp(ngvc,vgqc,ngvc,sfacgq)
78! generate the regularised Coulomb Green's function in G+q-space
79 call gengclgq(.true.,iq,ngvc,gqc,gclgq)
80! compute the required spherical Bessel functions
81 call genjlgprmt(lnpsd,ngvc,gqc,ngvc,jlgqrmt)
82! find the matching coefficients
83 call match(ngk(1,ik),vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
84! get the eigenvectors from file for non-reduced k-points
85 call getevecfv(filext,0,vkl(:,ik),vgkl(:,:,1,ik),evecfv)
86 call getevecsv(filext,0,vkl(:,ik),evecsv)
87! count and index occupied states
88 nst=0
89 do ist3=1,nstsv
90 if (evalsv(ist3,jk) > efermi) cycle
91 nst=nst+1
92 idx(nst)=ist3
93 end do
94! calculate the wavefunctions for occupied states
95 call genwfsv_sp(.false.,.false.,nst,idx,ngdgc,igfc,ngk(1,ik),igkig(:,1,ik), &
96 apwalm,evecfv,evecsv,wfmt2,ngtc,wfir2)
97 do ist3=1,nst
98! compute the complex overlap densities for all valence-valence states
99 do ist1=1,nstsv
100 call gencrho(.true.,.true.,ngtc,wfmt2(:,:,:,ist3),wfir2(:,:,ist3), &
101 wfmt1(:,:,:,ist1),wfir1(:,:,ist1),crhomt1(:,:,ist1),crhoir1(:,ist1))
102 end do
103! compute the complex overlap densities for all valence-core states
104 jc=0
105 do is=1,nspecies
106 nrc=nrcmt(is)
107 nrci=nrcmti(is)
108 npc=npcmt(is)
109 do ia=1,natoms(is)
110 ias=idxas(ia,is)
111 do ist1=1,nstsp(is)
112 if (spcore(ist1,is)) then
113 do m1=-ksp(ist1,is),ksp(ist1,is)-1
114 jc=jc+1
115! generate the core wavefunction in spherical coordinates (pass in m-1/2)
116 call wavefcr(.false.,lradstp,is,ia,ist1,m1,npcmtmax,wfcr1)
117 if (spinpol) then
118 call crho2(npc,wfmt2(:,ias,1,ist3),wfmt2(:,ias,2,ist3),wfcr1, &
119 wfcr1(:,2),crhomt2(:,jc))
120 else
121 call crho1(npc,wfmt2(:,ias,1,ist3),wfcr1,crhomt2(:,jc))
122 end if
123! convert to spherical harmonics
124 call cfshtip(nrc,nrci,crhomt2(:,jc))
125 end do
126 end if
127 end do
128 end do
129 end do
130 do ist2=1,nstsv
131 if (evalsv(ist2,ikp) > efermi) then
132! calculate the Coulomb potential
134 crhomt1(:,:,ist2),cvclmt)
136 gclgq,ngvc,jlgqrmt,ylmgq,sfacgq,crhoir1(:,ist2),npcmtmax,cvclmt,cvclir)
137 cvclir(:)=cvclir(:)*cfrc(:)
138!----------------------------------------------!
139! valence-valence-valence contribution !
140!----------------------------------------------!
141 do ist1=1,nstsv
142 if (evalsv(ist1,ikp) < efermi) then
143 z1=zcfinp(crhomt1(:,:,ist1),crhoir1(:,ist1),cvclmt,cvclir)
144 vclvv(ist1,ist2)=vclvv(ist1,ist2)-wqptnr*z1
145 end if
146 end do
147!-------------------------------------------!
148! core-valence-valence contribution !
149!-------------------------------------------!
150 jc=0
151 do is=1,nspecies
152 nrc=nrcmt(is)
153 nrci=nrcmti(is)
154 do ia=1,natoms(is)
155 ias=idxas(ia,is)
156 ic=0
157 do ist1=1,nstsp(is)
158 if (spcore(ist1,is)) then
159 do m1=-ksp(ist1,is),ksp(ist1,is)-1
160 ic=ic+1
161 jc=jc+1
162 z1=zcfmtinp(nrc,nrci,wr2cmt(:,is),crhomt2(:,jc),cvclmt(:,ias))
163 vclcv(ic,ias,ist2)=vclcv(ic,ias,ist2)-wqptnr*z1
164 end do
165! end loop over ist1
166 end if
167 end do
168! end loops over atoms and species
169 end do
170 end do
171! end loop over ist2
172 end if
173 end do
174! end loop over ist3
175 end do
176! end loop over non-reduced k-point set
177end do
178! begin loops over atoms and species
179do is=1,nspecies
180 nrc=nrcmt(is)
181 nrci=nrcmti(is)
182 npc=npcmt(is)
183 do ia=1,natoms(is)
184 ias=idxas(ia,is)
185 do ist3=1,nstsp(is)
186 if (spcore(ist3,is)) then
187 do m1=-ksp(ist3,is),ksp(ist3,is)-1
188! generate the core wavefunction in spherical coordinates (pass in m-1/2)
189 call wavefcr(.false.,lradstp,is,ia,ist3,m1,npcmtmax,wfcr1)
190! compute the complex overlap densities for the core-valence states
191 do ist1=1,nstsv
192 if (spinpol) then
193 call crho2(npc,wfcr1,wfcr1(:,2),wfmt1(:,ias,1,ist1), &
194 wfmt1(:,ias,2,ist1),crhomt1(:,ias,ist1))
195 else
196 call crho1(npc,wfcr1,wfmt1(:,ias,1,ist1),crhomt1(:,ias,ist1))
197 end if
198 call cfshtip(nrc,nrci,crhomt1(:,ias,ist1))
199 end do
200! compute the complex overlap densities for the core-core states
201 ic=0
202 do ist1=1,nstsp(is)
203 if (spcore(ist1,is)) then
204 do m2=-ksp(ist1,is),ksp(ist1,is)-1
205 ic=ic+1
206 call wavefcr(.false.,lradstp,is,ia,ist1,m2,npcmtmax,wfcr2)
207 call crho2(npc,wfcr1,wfcr1(:,2),wfcr2,wfcr2(:,2),crhomt2(:,ic))
208 call cfshtip(nrc,nrci,crhomt2(:,ic))
209 end do
210 end if
211 end do
212 do ist2=1,nstsv
213 if (evalsv(ist2,ikp) > efermi) then
214! calculate the Coulomb potential
215 call cpotclmt(nrc,nrci,nrcmtmax,rlcmt(:,:,is),wprcmt(:,:,is), &
216 crhomt1(:,ias,ist2),cvclmt)
217!-------------------------------------------!
218! valence-core-valence contribution !
219!-------------------------------------------!
220 do ist1=1,nstsv
221 if (evalsv(ist1,ikp) < efermi) then
222 z1=zcfmtinp(nrc,nrci,wr2cmt(:,is),crhomt1(:,ias,ist1),cvclmt)
223 vclvv(ist1,ist2)=vclvv(ist1,ist2)-z1
224 end if
225 end do
226!----------------------------------------!
227! core-core-valence contribution !
228!----------------------------------------!
229 ic=0
230 do ist1=1,nstsp(is)
231 if (spcore(ist1,is)) then
232 do m2=-ksp(ist1,is),ksp(ist1,is)-1
233 ic=ic+1
234 z1=zcfmtinp(nrc,nrci,wr2cmt(:,is),crhomt2(:,ic),cvclmt)
235 vclcv(ic,ias,ist2)=vclcv(ic,ias,ist2)-z1
236 end do
237! end loop over ist1
238 end if
239 end do
240! end loop over ist2
241 end if
242 end do
243! end loops over ist3 and m1
244 end do
245 end if
246 end do
247! end loops over atoms and species
248 end do
249end do
250deallocate(vgqc,gqc,gclgq,jlgqrmt)
251deallocate(apwalm,evecfv,evecsv,ylmgq,sfacgq)
252deallocate(wfmt1,wfir1,wfmt2,wfir2,wfcr1,wfcr2)
253deallocate(crhomt1,crhomt2,crhoir1)
254deallocate(cvclmt,cvclir)
255return
256
257contains
258
259pure subroutine crho1(n,wf1,wf2,crho)
260implicit none
261integer, intent(in) :: n
262complex(4), intent(in) :: wf1(n),wf2(n)
263complex(4), intent(out) :: crho(n)
264crho(1:n)=conjg(wf1(1:n))*wf2(1:n)
265end subroutine
266
267pure subroutine crho2(n,wf11,wf12,wf21,wf22,crho)
268implicit none
269integer, intent(in) :: n
270complex(4), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
271complex(4), intent(out) :: crho(n)
272crho(1:n)=conjg(wf11(1:n))*wf21(1:n)+conjg(wf12(1:n))*wf22(1:n)
273end subroutine
274
275end subroutine
276!EOC
277
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 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 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 nstcr
Definition modmain.f90:129
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
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
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), 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
subroutine oepvclk(ikp, vclcv, vclvv)
Definition oepvclk.f90:7
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