The Elk Code
 
Loading...
Searching...
No Matches
genspfxcg.f90
Go to the documentation of this file.
1
2! Copyright (C) 2013 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 genspfxcg(fxc)
7use modmain
8implicit none
9! arguments
10complex(8), intent(out) :: fxc(ngrf,4,ngrf,4)
11! local variables
12integer ig,jg,kg
13integer iv(3),i,j
14complex(8) z1
15! allocatable arrays
16real(8), allocatable :: fxcmt(:,:,:,:),fxcir(:,:,:)
17complex(8), allocatable :: fxcg(:)
18allocate(fxcmt(npmtmax,natmtot,4,4),fxcir(ngtot,4,4))
19allocate(fxcg(ngvec))
20! generate the kernel f_xc in real-space
21call genspfxcr(.true.,fxcmt,fxcir)
22! Fourier transform the kernel to G-space
23do i=1,4
24 do j=i,4
25 call zftrf(ngvec,ivg,vgc,fxcmt(:,:,i,j),fxcir(:,i,j),fxcg)
26 do ig=1,ngrf
27 do jg=1,ngrf
28 iv(1:3)=ivg(1:3,ig)-ivg(1:3,jg)
29 if ((iv(1) >= intgv(1,1)).and.(iv(1) <= intgv(2,1)).and. &
30 (iv(2) >= intgv(1,2)).and.(iv(2) <= intgv(2,2)).and. &
31 (iv(3) >= intgv(1,3)).and.(iv(3) <= intgv(2,3))) then
32 kg=ivgig(iv(1),iv(2),iv(3))
33 if (kg > ngvec) cycle
34 z1=fxcg(kg)
35 fxc(ig,i,jg,j)=z1
36 fxc(jg,j,ig,i)=conjg(z1)
37 end if
38 end do
39 end do
40 end do
41end do
42deallocate(fxcmt,fxcir,fxcg)
43end subroutine
44
subroutine genspfxcg(fxc)
Definition genspfxcg.f90:7
subroutine genspfxcr(tsh, fxcmt, fxcir)
Definition genspfxcr.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