The Elk Code
 
Loading...
Searching...
No Matches
ggair_sp_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_sp_2b
8! !INTERFACE:
9subroutine ggair_sp_2b(g2up,g2dn,gvup,gvdn,vxup,vxdn,vcup,vcdn,dxdgu2,dxdgd2, &
10 dxdgud,dcdgu2,dcdgd2,dcdgud)
11! !USES:
12use modmain
13! !DESCRIPTION:
14! Post processing step of interstitial gradients for GGA type 2. See routine
15! {\tt ggamt\_sp\_2a} for full details.
16!
17! !REVISION HISTORY:
18! Created November 2009 (JKD and TMcQ)
19!EOP
20!BOC
21implicit none
22real(8), intent(in) :: g2up(ngtot),g2dn(ngtot)
23real(8), intent(in) :: gvup(ngtot,3),gvdn(ngtot,3)
24real(8), intent(inout) :: vxup(ngtot),vxdn(ngtot)
25real(8), intent(inout) :: vcup(ngtot),vcdn(ngtot)
26real(8), intent(in) :: dxdgu2(ngtot),dxdgd2(ngtot),dxdgud(ngtot)
27real(8), intent(in) :: dcdgu2(ngtot),dcdgd2(ngtot),dcdgud(ngtot)
28! local variables
29integer ig,ifg,i
30! allocatable arrays
31real(8), allocatable :: rfir1(:),rfir2(:)
32complex(8), allocatable :: zfft1(:),zfft2(:)
33allocate(rfir1(ngtot),rfir2(ngtot))
34allocate(zfft1(nfgrz),zfft2(nfgrz))
35!------------------!
36! exchange !
37!------------------!
38! compute grad dxdgu2
39call rzfftifc(3,ngridg,-1,dxdgu2,zfft1)
40! (grad dxdgu2).(grad rhoup)
41rfir1(:)=0.d0
42do i=1,3
43 do ifg=1,nfgrz
44 ig=igrzf(ifg)
45 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
46 end do
47 call rzfftifc(3,ngridg,1,rfir2,zfft2)
48 rfir1(:)=rfir1(:)+rfir2(:)*gvup(:,i)
49end do
50vxup(:)=vxup(:)-2.d0*(rfir1(:)+dxdgu2(:)*g2up(:))-dxdgud(:)*g2dn(:)
51! compute grad dxdgd2
52call rzfftifc(3,ngridg,-1,dxdgd2,zfft1)
53! (grad dxdgd2).(grad rhodn)
54rfir1(:)=0.d0
55do i=1,3
56 do ifg=1,nfgrz
57 ig=igrzf(ifg)
58 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
59 end do
60 call rzfftifc(3,ngridg,1,rfir2,zfft2)
61 rfir1(:)=rfir1(:)+rfir2(:)*gvdn(:,i)
62end do
63vxdn(:)=vxdn(:)-2.d0*(rfir1(:)+dxdgd2(:)*g2dn(:))-dxdgud(:)*g2up(:)
64! compute grad dxdgud
65call rzfftifc(3,ngridg,-1,dxdgud,zfft1)
66! (grad dxdgud).(grad rhodn) and (grad dxdgud).(grad rhoup)
67do i=1,3
68 do ifg=1,nfgrz
69 ig=igrzf(ifg)
70 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
71 end do
72 call rzfftifc(3,ngridg,1,rfir2,zfft2)
73 vxup(:)=vxup(:)-rfir2(:)*gvdn(:,i)
74 vxdn(:)=vxdn(:)-rfir2(:)*gvup(:,i)
75end do
76!---------------------!
77! correlation !
78!---------------------!
79! compute grad dcdgu2
80call rzfftifc(3,ngridg,-1,dcdgu2,zfft1)
81! (grad dcdgu2).(grad rhoup)
82rfir1(:)=0.d0
83do i=1,3
84 do ifg=1,nfgrz
85 ig=igrzf(ifg)
86 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
87 end do
88 call rzfftifc(3,ngridg,1,rfir2,zfft2)
89 rfir1(:)=rfir1(:)+rfir2(:)*gvup(:,i)
90end do
91vcup(:)=vcup(:)-2.d0*(rfir1(:)+dcdgu2(:)*g2up(:))-dcdgud(:)*g2dn(:)
92! compute grad dcdgd2
93call rzfftifc(3,ngridg,-1,dcdgd2,zfft1)
94! (grad dcdgd2).(grad rhodn)
95rfir1(:)=0.d0
96do i=1,3
97 do ifg=1,nfgrz
98 ig=igrzf(ifg)
99 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
100 end do
101 call rzfftifc(3,ngridg,1,rfir2,zfft2)
102 rfir1(:)=rfir1(:)+rfir2(:)*gvdn(:,i)
103end do
104vcdn(:)=vcdn(:)-2.d0*(rfir1(:)+dcdgd2(:)*g2dn(:))-dcdgud(:)*g2up(:)
105! compute grad dcdgud
106call rzfftifc(3,ngridg,-1,dcdgud,zfft1)
107! (grad dcdgud).(grad rhodn) and (grad dcdgud).(grad rhoup)
108do i=1,3
109 do ifg=1,nfgrz
110 ig=igrzf(ifg)
111 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
112 end do
113 call rzfftifc(3,ngridg,1,rfir2,zfft2)
114 vcup(:)=vcup(:)-rfir2(:)*gvdn(:,i)
115 vcdn(:)=vcdn(:)-rfir2(:)*gvup(:,i)
116end do
117deallocate(rfir1,rfir2,zfft1,zfft2)
118end subroutine
119!EOC
120
subroutine ggair_sp_2b(g2up, g2dn, gvup, gvdn, vxup, vxdn, vcup, vcdn, dxdgu2, dxdgd2, dxdgud, dcdgu2, dcdgd2, dcdgud)
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)