The Elk Code
ggair_sp_1.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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_1
8 ! !INTERFACE:
9 subroutine ggair_sp_1(rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! rhoup : spin-up density (in,real(ngtot))
12 ! rhodn : spin-down density (in,real(ngtot))
13 ! grho : |grad rho| (out,real(ngtot))
14 ! gup : |grad rhoup| (out,real(ngtot))
15 ! gdn : |grad rhodn| (out,real(ngtot))
16 ! g2up : grad^2 rhoup (out,real(ngtot))
17 ! g2dn : grad^2 rhodn (out,real(ngtot))
18 ! g3rho : (grad rho).(grad |grad rho|) (out,real(ngtot))
19 ! g3up : (grad rhoup).(grad |grad rhoup|) (out,real(ngtot))
20 ! g3dn : (grad rhodn).(grad |grad rhodn|) (out,real(ngtot))
21 ! !DESCRIPTION:
22 ! Computes $|\nabla\rho|$, $|\nabla\rho^{\uparrow}|$,
23 ! $|\nabla\rho^{\downarrow}|$, $\nabla^2\rho^{\uparrow}$,
24 ! $\nabla^2\rho^{\downarrow}$, $\nabla\rho\cdot(\nabla|\nabla\rho|)$,
25 ! $\nabla\rho^{\uparrow}\cdot(\nabla|\nabla\rho^{\uparrow}|)$ and
26 ! $\nabla\rho^{\downarrow}\cdot(\nabla|\nabla\rho^{\downarrow}|)$ for the
27 ! interstitial charge density, as required by the generalised gradient
28 ! approximation functionals of type 1 for spin-polarised densities. See
29 ! routines {\tt potxc} and {\tt modxcifc}.
30 !
31 ! !REVISION HISTORY:
32 ! Created October 2004 (JKD)
33 ! Simplified and improved, October 2009 (JKD)
34 !EOP
35 !BOC
36 use modmain
37 implicit none
38 ! arguments
39 real(8), intent(in) :: rhoup(ngtot),rhodn(ngtot)
40 real(8), intent(out) :: grho(ngtot),gup(ngtot),gdn(ngtot)
41 real(8), intent(out) :: g2up(ngtot),g2dn(ngtot)
42 real(8), intent(out) :: g3rho(ngtot),g3up(ngtot),g3dn(ngtot)
43 ! local variables
44 integer ig,ifg,i
45 ! allocatable arrays
46 real(8), allocatable :: gvup(:,:),gvdn(:,:),rfir(:)
47 complex(8), allocatable :: zfft1(:),zfft2(:)
48 allocate(gvup(ngtot,3),gvdn(ngtot,3),rfir(ngtot))
49 allocate(zfft1(nfgrz),zfft2(nfgrz))
50 !------------!
51 ! ρ↑ !
52 !------------!
53 call rzfftifc(3,ngridg,-1,rhoup,zfft1)
54 ! |∇ρ↑|
55 do i=1,3
56  do ifg=1,nfgrz
57  ig=igrzf(ifg)
58  if (ig <= ngvc) then
59  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
60  else
61  zfft2(ifg)=0.d0
62  end if
63  end do
64  call rzfftifc(3,ngridg,1,gvup(:,i),zfft2)
65 end do
66 gup(:)=sqrt(gvup(:,1)**2+gvup(:,2)**2+gvup(:,3)**2)
67 ! ∇²ρ↑
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,g2up,zfft2)
77 ! (∇ρ↑)⋅(∇|∇ρ↑|)
78 call rzfftifc(3,ngridg,-1,gup,zfft1)
79 g3up(:)=0.d0
80 do i=1,3
81  do ifg=1,nfgrz
82  ig=igrzf(ifg)
83  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
84  end do
85  call rzfftifc(3,ngridg,1,rfir,zfft2)
86  g3up(:)=g3up(:)+gvup(:,i)*rfir(:)
87 end do
88 !------------!
89 ! ρ↓ !
90 !------------!
91 call rzfftifc(3,ngridg,-1,rhodn,zfft1)
92 ! |∇ρ↓|
93 do i=1,3
94  do ifg=1,nfgrz
95  ig=igrzf(ifg)
96  if (ig <= ngvc) then
97  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
98  else
99  zfft2(ifg)=0.d0
100  end if
101  end do
102  call rzfftifc(3,ngridg,1,gvdn(:,i),zfft2)
103 end do
104 gdn(:)=sqrt(gvdn(:,1)**2+gvdn(:,2)**2+gvdn(:,3)**2)
105 ! ∇²ρ↓
106 do ifg=1,nfgrz
107  ig=igrzf(ifg)
108  if (ig <= ngvc) then
109  zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
110  else
111  zfft2(ifg)=0.d0
112  end if
113 end do
114 call rzfftifc(3,ngridg,1,g2dn,zfft2)
115 ! (∇ρ↓)⋅(∇|∇ρ↓|)
116 call rzfftifc(3,ngridg,-1,gdn,zfft1)
117 g3dn(:)=0.d0
118 do i=1,3
119  do ifg=1,nfgrz
120  ig=igrzf(ifg)
121  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
122  end do
123  call rzfftifc(3,ngridg,1,rfir,zfft2)
124  g3dn(:)=g3dn(:)+gvdn(:,i)*rfir(:)
125 end do
126 !-----------!
127 ! ρ !
128 !-----------!
129 ! |∇ρ|
130 grho(:)=sqrt((gvup(:,1)+gvdn(:,1))**2 &
131  +(gvup(:,2)+gvdn(:,2))**2 &
132  +(gvup(:,3)+gvdn(:,3))**2)
133 ! (∇ρ)⋅(∇|∇ρ|)
134 call rzfftifc(3,ngridg,-1,grho,zfft1)
135 g3rho(:)=0.d0
136 do i=1,3
137  do ifg=1,nfgrz
138  ig=igrzf(ifg)
139  zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
140  end do
141  call rzfftifc(3,ngridg,1,rfir,zfft2)
142  g3rho(:)=g3rho(:)+(gvup(:,i)+gvdn(:,i))*rfir(:)
143 end do
144 deallocate(gvup,gvdn,rfir,zfft1,zfft2)
145 end subroutine
146 !EOC
147 
integer, dimension(3) ngridg
Definition: modmain.f90:386
subroutine ggair_sp_1(rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn)
Definition: ggair_sp_1.f90:10
integer nfgrz
Definition: modmain.f90:412
integer ngvc
Definition: modmain.f90:398
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer, dimension(:), allocatable igrzf
Definition: modmain.f90:416
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
subroutine rzfftifc(nd, n, sgn, r, z)