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