The Elk Code
Loading...
Searching...
No Matches
genvfxcg.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2012 S. Sharma, J. K. Dewhurst 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
genvfxcg
(gclgq,nm,vfxc)
7
use
modmain
8
implicit none
9
! arguments
10
real
(8),
intent(in)
:: gclgq(ngrf)
11
integer
,
intent(in)
:: nm
12
complex(8)
,
intent(out)
:: vfxc(nm,nm,nwrf)
13
! local variables
14
integer
ig,jg,kg,i1,i2,i3
15
complex(8)
z1
16
! allocatable arrays
17
real
(8),
allocatable
:: fxcmt(:,:),fxcir(:)
18
complex(8)
,
allocatable
:: fxcg(:)
19
allocate
(fxcmt(
npmtmax
,
natmtot
),fxcir(
ngtot
))
20
allocate
(fxcg(
ngvec
))
21
! generate the kernel f_xc in real-space
22
call
genfxcr
(.true.,fxcmt,fxcir)
23
! Fourier transform the kernel to G-space
24
call
zftrf
(
ngvec
,
ivg
,
vgc
,fxcmt,fxcir,fxcg)
25
do
ig=1,ngrf
26
do
jg=1,ngrf
27
i1=
ivg
(1,ig)-
ivg
(1,jg)
28
i2=
ivg
(2,ig)-
ivg
(2,jg)
29
i3=
ivg
(3,ig)-
ivg
(3,jg)
30
if
((i1 >=
intgv
(1,1)).and.(i1 <=
intgv
(2,1)).and. &
31
(i2 >=
intgv
(1,2)).and.(i2 <=
intgv
(2,2)).and. &
32
(i3 >=
intgv
(1,3)).and.(i3 <=
intgv
(2,3)))
then
33
kg=
ivgig
(i1,i2,i3)
34
if
(kg <=
ngvec
)
then
35
z1=fxcg(kg)/(gclgq(ig)*gclgq(jg))
36
vfxc(ig,jg,1:nwrf)=z1
37
else
38
vfxc(ig,jg,1:nwrf)=0.d0
39
end if
40
end if
41
end do
42
end do
43
deallocate
(fxcmt,fxcir,fxcg)
44
end subroutine
45
genfxcr
subroutine genfxcr(tsh, fxcmt, fxcir)
Definition
genfxcr.f90:7
genvfxcg
subroutine genvfxcg(gclgq, nm, vfxc)
Definition
genvfxcg.f90:7
modmain
Definition
modmain.f90:6
modmain::intgv
integer, dimension(2, 3) intgv
Definition
modmain.f90:394
modmain::ngtot
integer ngtot
Definition
modmain.f90:390
modmain::ngvec
integer ngvec
Definition
modmain.f90:396
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::npmtmax
integer npmtmax
Definition
modmain.f90:216
modmain::ivg
integer, dimension(:,:), allocatable ivg
Definition
modmain.f90:400
modmain::vgc
real(8), dimension(:,:), allocatable vgc
Definition
modmain.f90:420
modmain::ivgig
integer, dimension(:,:,:), allocatable ivgig
Definition
modmain.f90:402
zftrf
subroutine zftrf(npv, ivp, vpc, rfmt, rfir, zfp)
Definition
zftrf.f90:10
genvfxcg.f90
Generated by
1.9.8