The Elk Code
Loading...
Searching...
No Matches
ggair_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_2a
8
! !INTERFACE:
9
subroutine
ggair_2a
(rho,g2rho,gvrho,grho2)
10
! !USES:
11
use
modmain
12
! !DESCRIPTION:
13
! Spin-unpolarised version of {\tt ggair\_sp\_2a}.
14
!
15
! !REVISION HISTORY:
16
! Created November 2009 (JKD and TMcQ)
17
!EOP
18
!BOC
19
implicit none
20
! arguments
21
real
(8),
intent(in)
:: rho(ngtot)
22
real
(8),
intent(out)
:: g2rho(ngtot),gvrho(ngtot,3),grho2(ngtot)
23
! local variables
24
integer
ig,ifg,i
25
! allocatable arrays
26
complex(8)
,
allocatable
:: zfft1(:),zfft2(:)
27
allocate
(zfft1(
nfgrz
),zfft2(
nfgrz
))
28
! Fourier transform density to G-space
29
call
rzfftifc
(3,
ngridg
,-1,rho,zfft1)
30
! grad^2 rho
31
do
ifg=1,
nfgrz
32
ig=
igrzf
(ifg)
33
if
(ig <=
ngvc
)
then
34
zfft2(ifg)=-(
gc
(ig)**2)*zfft1(ifg)
35
else
36
zfft2(ifg)=0.d0
37
end if
38
end do
39
call
rzfftifc
(3,
ngridg
,1,g2rho,zfft2)
40
! grad rho
41
do
i=1,3
42
do
ifg=1,
nfgrz
43
ig=
igrzf
(ifg)
44
if
(ig <=
ngvc
)
then
45
zfft2(ifg)=
vgc
(i,ig)*
zi
*zfft1(ifg)
46
else
47
zfft2(ifg)=0.d0
48
end if
49
end do
50
call
rzfftifc
(3,
ngridg
,1,gvrho(:,i),zfft2)
51
end do
52
! (grad rho)^2
53
grho2(:)=gvrho(:,1)**2+gvrho(:,2)**2+gvrho(:,3)**2
54
deallocate
(zfft1,zfft2)
55
end subroutine
56
!EOC
57
ggair_2a
subroutine ggair_2a(rho, g2rho, gvrho, grho2)
Definition
ggair_2a.f90:10
modmain
Definition
modmain.f90:6
modmain::ngridg
integer, dimension(3) ngridg
Definition
modmain.f90:386
modmain::nfgrz
integer nfgrz
Definition
modmain.f90:412
modmain::igrzf
integer, dimension(:), allocatable igrzf
Definition
modmain.f90:416
modmain::zi
complex(8), parameter zi
Definition
modmain.f90:1239
modmain::vgc
real(8), dimension(:,:), allocatable vgc
Definition
modmain.f90:420
modmain::ngvc
integer ngvc
Definition
modmain.f90:398
modmain::gc
real(8), dimension(:), allocatable gc
Definition
modmain.f90:422
rzfftifc
subroutine rzfftifc(nd, n, sgn, r, z)
Definition
zfftifc_fftw.f90:32
ggair_2a.f90
Generated by
1.9.8