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
6subroutine genvfxcg(gclgq,nm,vfxc)
7use modmain
8implicit none
9! arguments
10real(8), intent(in) :: gclgq(ngrf)
11integer, intent(in) :: nm
12complex(8), intent(out) :: vfxc(nm,nm,nwrf)
13! local variables
14integer ig,jg,kg,i1,i2,i3
15complex(8) z1
16! allocatable arrays
17real(8), allocatable :: fxcmt(:,:),fxcir(:)
18complex(8), allocatable :: fxcg(:)
19allocate(fxcmt(npmtmax,natmtot),fxcir(ngtot))
20allocate(fxcg(ngvec))
21! generate the kernel f_xc in real-space
22call genfxcr(.true.,fxcmt,fxcir)
23! Fourier transform the kernel to G-space
24call zftrf(ngvec,ivg,vgc,fxcmt,fxcir,fxcg)
25do 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
42end do
43deallocate(fxcmt,fxcir,fxcg)
44end subroutine
45
subroutine genfxcr(tsh, fxcmt, fxcir)
Definition genfxcr.f90:7
subroutine genvfxcg(gclgq, nm, vfxc)
Definition genvfxcg.f90:7
integer, dimension(2, 3) intgv
Definition modmain.f90:394
integer ngtot
Definition modmain.f90:390
integer ngvec
Definition modmain.f90:396
integer natmtot
Definition modmain.f90:40
integer npmtmax
Definition modmain.f90:216
integer, dimension(:,:), allocatable ivg
Definition modmain.f90:400
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer, dimension(:,:,:), allocatable ivgig
Definition modmain.f90:402
subroutine zftrf(npv, ivp, vpc, rfmt, rfir, zfp)
Definition zftrf.f90:10