The Elk Code
 
Loading...
Searching...
No Matches
ggair_2b.f90
Go to the documentation of this file.
1
2! Copyright (C) 2009 T. McQueen and J. K. Dewhurst.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: ggair_2b
8! !INTERFACE:
9subroutine ggair_2b(g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Spin-unpolarised version of {\tt ggair\_sp\_2b}.
14!
15! !REVISION HISTORY:
16! Created November 2009 (JKD and TMcQ)
17!EOP
18!BOC
19implicit none
20! arguments
21real(8), intent(in) :: g2rho(ngtot),gvrho(ngtot,3)
22real(8), intent(inout) :: vx(ngtot),vc(ngtot)
23real(8), intent(in) :: dxdgr2(ngtot),dcdgr2(ngtot)
24! local variables
25integer ig,ifg,i
26! allocatable arrays
27real(8), allocatable :: rfir1(:),rfir2(:)
28complex(8), allocatable :: zfft1(:),zfft2(:)
29allocate(rfir1(ngtot),rfir2(ngtot))
30allocate(zfft1(nfgrz),zfft2(nfgrz))
31!------------------!
32! exchange !
33!------------------!
34! compute grad dxdgr2
35call rzfftifc(3,ngridg,-1,dxdgr2,zfft1)
36! (grad dxdgr2).(grad rho)
37rfir1(:)=0.d0
38do i=1,3
39 do ifg=1,nfgrz
40 ig=igrzf(ifg)
41 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
42 end do
43 call rzfftifc(3,ngridg,1,rfir2,zfft2)
44 rfir1(:)=rfir1(:)+rfir2(:)*gvrho(:,i)
45end do
46vx(:)=vx(:)-2.d0*(rfir1(:)+dxdgr2(:)*g2rho(:))
47!---------------------!
48! correlation !
49!---------------------!
50! compute grad dcdgr2
51call rzfftifc(3,ngridg,-1,dcdgr2,zfft1)
52! (grad dcdgr2).(grad rho)
53rfir1(:)=0.d0
54do i=1,3
55 do ifg=1,nfgrz
56 ig=igrzf(ifg)
57 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
58 end do
59 call rzfftifc(3,ngridg,1,rfir2,zfft2)
60 rfir1(:)=rfir1(:)+rfir2(:)*gvrho(:,i)
61end do
62vc(:)=vc(:)-2.d0*(rfir1(:)+dcdgr2(:)*g2rho(:))
63deallocate(rfir1,rfir2,zfft1,zfft2)
64end subroutine
65!EOC
66
subroutine ggair_2b(g2rho, gvrho, vx, vc, dxdgr2, dcdgr2)
Definition ggair_2b.f90:10
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)