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 :: rfir(:) 32complex(8), allocatable :: zfft1(:),zfft2(:) 33allocate(rfir(ngtot)) 34allocate(zfft1(ngtot),zfft2(ngtot)) 35!------------------! 36! exchange ! 37!------------------! 38! compute grad dxdgu2 39zfft1(:)=dxdgu2(:) 40call zfftifc(3,ngridg,-1,zfft1) 41! (grad dxdgu2).(grad rhoup) 42rfir(:)=0.d0 43do i=1,3 44 zfft2(:)=0.d0 45 do ig=1,ngvec 46 ifg=igfft(ig) 47 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 48 end do 49 call zfftifc(3,ngridg,1,zfft2) 50 rfir(:)=rfir(:)+dble(zfft2(:))*gvup(:,i) 51end do 52vxup(:)=vxup(:)-2.d0*(rfir(:)+dxdgu2(:)*g2up(:))-dxdgud(:)*g2dn(:) 53! compute grad dxdgd2 54zfft1(:)=dxdgd2(:) 55call zfftifc(3,ngridg,-1,zfft1) 56! (grad dxdgd2).(grad rhodn) 57rfir(:)=0.d0 58do i=1,3 59 zfft2(:)=0.d0 60 do ig=1,ngvec 61 ifg=igfft(ig) 62 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 63 end do 64 call zfftifc(3,ngridg,1,zfft2) 65 rfir(:)=rfir(:)+dble(zfft2(:))*gvdn(:,i) 66end do 67vxdn(:)=vxdn(:)-2.d0*(rfir(:)+dxdgd2(:)*g2dn(:))-dxdgud(:)*g2up(:) 68! compute grad dxdgud 69zfft1(:)=dxdgud(:) 70call zfftifc(3,ngridg,-1,zfft1) 71! (grad dxdgud).(grad rhodn) and (grad dxdgud).(grad rhoup) 72do i=1,3 73 zfft2(:)=0.d0 74 do ig=1,ngvec 75 ifg=igfft(ig) 76 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 77 end do 78 call zfftifc(3,ngridg,1,zfft2) 79 vxup(:)=vxup(:)-dble(zfft2(:))*gvdn(:,i) 80 vxdn(:)=vxdn(:)-dble(zfft2(:))*gvup(:,i) 81end do 82!---------------------! 83! correlation ! 84!---------------------! 85! compute grad dcdgu2 86zfft1(:)=dcdgu2(:) 87call zfftifc(3,ngridg,-1,zfft1) 88! (grad dcdgu2).(grad rhoup) 89rfir(:)=0.d0 90do i=1,3 91 zfft2(:)=0.d0 92 do ig=1,ngvec 93 ifg=igfft(ig) 94 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 95 end do 96 call zfftifc(3,ngridg,1,zfft2) 97 rfir(:)=rfir(:)+dble(zfft2(:))*gvup(:,i) 98end do 99vcup(:)=vcup(:)-2.d0*(rfir(:)+dcdgu2(:)*g2up(:))-dcdgud(:)*g2dn(:) 100! compute grad dcdgd2 101zfft1(:)=dcdgd2(:) 102call zfftifc(3,ngridg,-1,zfft1) 103! (grad dcdgd2).(grad rhodn) 104rfir(:)=0.d0 105do i=1,3 106 zfft2(:)=0.d0 107 do ig=1,ngvec 108 ifg=igfft(ig) 109 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 110 end do 111 call zfftifc(3,ngridg,1,zfft2) 112 rfir(:)=rfir(:)+dble(zfft2(:))*gvdn(:,i) 113end do 114vcdn(:)=vcdn(:)-2.d0*(rfir(:)+dcdgd2(:)*g2dn(:))-dcdgud(:)*g2up(:) 115! compute grad dcdgud 116zfft1(:)=dcdgud(:) 117call zfftifc(3,ngridg,-1,zfft1) 118! (grad dcdgud).(grad rhodn) and (grad dcdgud).(grad rhoup) 119do i=1,3 120 zfft2(:)=0.d0 121 do ig=1,ngvec 122 ifg=igfft(ig) 123 zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8) 124 end do 125 call zfftifc(3,ngridg,1,zfft2) 126 vcup(:)=vcup(:)-dble(zfft2(:))*gvdn(:,i) 127 vcdn(:)=vcdn(:)-dble(zfft2(:))*gvup(:,i) 128end do 129deallocate(rfir,zfft1,zfft2) 130end subroutine 131!EOC 132 133