The Elk Code
ggair_sp_2a.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_2a
8 ! !INTERFACE:
9 subroutine ggair_sp_2a(rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Computes the interstitial gradients $\nabla^2\rho^{\uparrow}$,
14 ! $\nabla^2\rho^{\downarrow}$, $\nabla\rho^{\uparrow}$,
15 ! $\nabla\rho^{\downarrow}$, $(\nabla\rho^{\uparrow})^2$,
16 ! $(\nabla\rho^{\downarrow})^2$ and
17 ! $\nabla\rho^{\uparrow}\cdot\nabla\rho^{\downarrow}$. These are used for GGA
18 ! functionals of type 2 and meta-GGA. See {\tt ggamt\_sp\_2a} for details.
19 !
20 ! !REVISION HISTORY:
21 ! Created November 2009 (JKD and TMcQ)
22 !EOP
23 !BOC
24 implicit none
25 ! arguments
26 real(8), intent(in) :: rhoup(ngtot),rhodn(ngtot)
27 real(8), intent(out) :: g2up(ngtot),g2dn(ngtot)
28 real(8), intent(out) :: gvup(ngtot,3),gvdn(ngtot,3)
29 real(8), intent(out) :: gup2(ngtot),gdn2(ngtot),gupdn(ngtot)
30 ! local variables
31 integer ig,ifg,i
32 ! allocatable arrays
33 complex(8), allocatable :: zfft1(:),zfft2(:)
34 allocate(zfft1(nfgrz),zfft2(nfgrz))
35 !------------!
36 ! ρ↑ !
37 !------------!
38 call rzfftifc(3,ngridg,-1,rhoup,zfft1)
39 ! compute ∇²ρ↑
40 do ifg=1,nfgrz
41  ig=igrzf(ifg)
42  if (ig <= ngvc) then
43  zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
44  else
45  zfft2(ifg)=0.d0
46  end if
47 end do
48 call rzfftifc(3,ngridg,1,g2up,zfft2)
49 ! ∇ρ↑
50 do i=1,3
51  do ifg=1,nfgrz
52  ig=igrzf(ifg)
53  if (ig <= ngvc) then
54  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
55  else
56  zfft2(ifg)=0.d0
57  end if
58  end do
59  call rzfftifc(3,ngridg,1,gvup(:,i),zfft2)
60 end do
61 ! (∇ρ↑)²
62 gup2(:)=gvup(:,1)**2+gvup(:,2)**2+gvup(:,3)**2
63 !------------!
64 ! ρ↓ !
65 !------------!
66 call rzfftifc(3,ngridg,-1,rhodn,zfft1)
67 ! compute ∇²ρ↓
68 do ifg=1,nfgrz
69  ig=igrzf(ifg)
70  if (ig <= ngvc) then
71  zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
72  else
73  zfft2(ifg)=0.d0
74  end if
75 end do
76 call rzfftifc(3,ngridg,1,g2dn,zfft2)
77 ! ∇ρ↓
78 do i=1,3
79  do ifg=1,nfgrz
80  ig=igrzf(ifg)
81  if (ig <= ngvc) then
82  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
83  else
84  zfft2(ifg)=0.d0
85  end if
86  end do
87  call rzfftifc(3,ngridg,1,gvdn(:,i),zfft2)
88 end do
89 ! (∇ρ↓)²
90 gdn2(:)=gvdn(:,1)**2+gvdn(:,2)**2+gvdn(:,3)**2
91 ! (∇ρ↑)⋅(∇ρ↓)
92 gupdn(:)=gvup(:,1)*gvdn(:,1)+gvup(:,2)*gvdn(:,2)+gvup(:,3)*gvdn(:,3)
93 deallocate(zfft1,zfft2)
94 end subroutine
95 !EOC
96 
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer nfgrz
Definition: modmain.f90:412
integer ngvc
Definition: modmain.f90:398
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
subroutine ggair_sp_2a(rhoup, rhodn, g2up, g2dn, gvup, gvdn, gup2, gdn2, gupdn)
Definition: ggair_sp_2a.f90:10
integer, dimension(:), allocatable igrzf
Definition: modmain.f90:416
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
subroutine rzfftifc(nd, n, sgn, r, z)