The Elk Code
 
Loading...
Searching...
No Matches
ggair_4.f90
Go to the documentation of this file.
1
2! Copyright (C) 2021 J. K. Dewhurst, S. Sharma 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 ggair_4(gvrho,vx,vc,wx,wc,dtdr,dtdgr2)
7use modmain
8implicit none
9! arguments
10real(8), intent(in) :: gvrho(ngtot,3)
11real(8), intent(inout) :: vx(ngtot),vc(ngtot)
12real(8), intent(in) :: wx(ngtot),wc(ngtot)
13real(8), intent(in) :: dtdr(ngtot),dtdgr2(ngtot)
14! local variables
15integer ig,ifg,i
16! allocatable arrays
17real(8), allocatable :: rfir1(:),rfir2(:)
18complex(8), allocatable :: zfft(:)
19allocate(rfir1(ngtot),rfir2(ngtot),zfft(nfgrz))
20!------------------!
21! exchange !
22!------------------!
23vx(:)=vx(:)+wx(:)*dtdr(:)
24rfir1(:)=wx(:)*dtdgr2(:)
25do i=1,3
26 rfir2(:)=rfir1(:)*gvrho(:,i)
27 call rzfftifc(3,ngridg,-1,rfir2,zfft)
28 do ifg=1,nfgrz
29 ig=igrzf(ifg)
30 zfft(ifg)=vgc(i,ig)*zi*zfft(ifg)
31 end do
32 call rzfftifc(3,ngridg,1,rfir2,zfft)
33 vx(:)=vx(:)-2.d0*rfir2(:)
34end do
35!---------------------!
36! correlation !
37!---------------------!
38vc(:)=vc(:)+wc(:)*dtdr(:)
39rfir1(:)=wc(:)*dtdgr2(:)
40do i=1,3
41 rfir2(:)=rfir1(:)*gvrho(:,i)
42 call rzfftifc(3,ngridg,-1,rfir2,zfft)
43 do ifg=1,nfgrz
44 ig=igrzf(ifg)
45 zfft(ifg)=vgc(i,ig)*zi*zfft(ifg)
46 end do
47 call rzfftifc(3,ngridg,1,rfir2,zfft)
48 vc(:)=vc(:)-2.d0*rfir2(:)
49end do
50deallocate(rfir1,rfir2,zfft)
51end subroutine
52
subroutine ggair_4(gvrho, vx, vc, wx, wc, dtdr, dtdgr2)
Definition ggair_4.f90:7
integer, dimension(3) ngridg
Definition modmain.f90:386
integer nfgrz
Definition modmain.f90:412
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
complex(8), parameter zi
Definition modmain.f90:1239
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
subroutine rzfftifc(nd, n, sgn, r, z)