The Elk Code
ggair_3.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2023 J. K. Dewhurst and S. Sharma.
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 ggair_3(vx,vc,wx,wc,dtdg2r)
7 use modmain
8 implicit none
9 ! arguments
10 real(8), intent(inout) :: vx(ngtot),vc(ngtot)
11 real(8), intent(in) :: wx(ngtot),wc(ngtot)
12 real(8), intent(in) :: dtdg2r(ngtot)
13 ! local variables
14 integer ig,ifg
15 ! allocatable arrays
16 real(8), allocatable :: rfir(:)
17 complex(8), allocatable :: zfft(:)
18 allocate(rfir(ngtot),zfft(nfgrz))
19 !------------------!
20 ! exchange !
21 !------------------!
22 ! Fourier transform (wx * dtdg2r) to G-space
23 rfir(:)=wx(:)*dtdg2r(:)
24 call rzfftifc(3,ngridg,-1,rfir,zfft)
25 ! ∇² (wx * dtdg2r)
26 do ifg=1,nfgrz
27  ig=igrzf(ifg)
28  zfft(ifg)=-(gc(ig)**2)*zfft(ifg)
29 end do
30 call rzfftifc(3,ngridg,1,rfir,zfft)
31 vx(:)=vx(:)+rfir(:)
32 !---------------------!
33 ! correlation !
34 !---------------------!
35 ! Fourier transform (wc * dtdg2r) to G-space
36 rfir(:)=wc(:)*dtdg2r(:)
37 call rzfftifc(3,ngridg,-1,rfir,zfft)
38 ! ∇² (wc * dtdg2r)
39 do ifg=1,nfgrz
40  ig=igrzf(ifg)
41  zfft(ifg)=-(gc(ig)**2)*zfft(ifg)
42 end do
43 call rzfftifc(3,ngridg,1,rfir,zfft)
44 vc(:)=vc(:)+rfir(:)
45 deallocate(rfir,zfft)
46 end subroutine
47 
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer nfgrz
Definition: modmain.f90:412
subroutine ggair_3(vx, vc, wx, wc, dtdg2r)
Definition: ggair_3.f90:7
integer, dimension(:), allocatable igrzf
Definition: modmain.f90:416
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
subroutine rzfftifc(nd, n, sgn, r, z)