The Elk Code
 
Loading...
Searching...
No Matches
ggair_1.f90
Go to the documentation of this file.
1
2! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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_1
8! !INTERFACE:
9subroutine ggair_1(rho,grho,g2rho,g3rho)
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Spin-unpolarised version of {\tt ggair\_sp\_1}.
14!
15! !REVISION HISTORY:
16! Created November 2009 (JKD)
17!EOP
18!BOC
19implicit none
20! arguments
21real(8), intent(in) :: rho(ngtot)
22real(8), intent(out) :: grho(ngtot),g2rho(ngtot),g3rho(ngtot)
23! local variables
24integer ig,ifg,i
25! allocatable arrays
26real(8), allocatable :: gvrho(:,:),rfir(:)
27complex(8), allocatable :: zfft1(:),zfft2(:)
28allocate(gvrho(ngtot,3),rfir(ngtot))
29allocate(zfft1(nfgrz),zfft2(nfgrz))
30call rzfftifc(3,ngridg,-1,rho,zfft1)
31! |grad rho|
32do i=1,3
33 do ifg=1,nfgrz
34 ig=igrzf(ifg)
35 if (ig <= ngvc) then
36 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
37 else
38 zfft2(ifg)=0.d0
39 end if
40 end do
41 call rzfftifc(3,ngridg,1,gvrho(:,i),zfft2)
42end do
43grho(:)=sqrt(gvrho(:,1)**2+gvrho(:,2)**2+gvrho(:,3)**2)
44! grad^2 rho
45do ifg=1,nfgrz
46 ig=igrzf(ifg)
47 if (ig <= ngvc) then
48 zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
49 else
50 zfft2(ifg)=0.d0
51 end if
52end do
53call rzfftifc(3,ngridg,1,g2rho,zfft2)
54! (grad rho).(grad |grad rho|)
55call rzfftifc(3,ngridg,-1,grho,zfft1)
56g3rho(:)=0.d0
57do i=1,3
58 do ifg=1,nfgrz
59 ig=igrzf(ifg)
60 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
61 end do
62 call rzfftifc(3,ngridg,1,rfir,zfft2)
63 g3rho(:)=g3rho(:)+gvrho(:,i)*rfir(:)
64end do
65deallocate(gvrho,rfir,zfft1,zfft2)
66end subroutine
67!EOC
68
subroutine ggair_1(rho, grho, g2rho, g3rho)
Definition ggair_1.f90:10
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
integer ngvc
Definition modmain.f90:398
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
subroutine rzfftifc(nd, n, sgn, r, z)