File: ggair_1.f90

package info (click to toggle)
elkcode 10.8.1-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 10,672 kB
  • sloc: f90: 52,747; perl: 964; makefile: 352; sh: 296; python: 105; ansic: 67
file content (68 lines) | stat: -rw-r--r-- 1,622 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68

! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: ggair_1
! !INTERFACE:
subroutine ggair_1(rho,grho,g2rho,g3rho)
! !USES:
use modmain
! !DESCRIPTION:
!   Spin-unpolarised version of {\tt ggair\_sp\_1}.
!
! !REVISION HISTORY:
!   Created November 2009 (JKD)
!EOP
!BOC
implicit none
! arguments
real(8), intent(in) :: rho(ngtot)
real(8), intent(out) :: grho(ngtot),g2rho(ngtot),g3rho(ngtot)
! local variables
integer ig,ifg,i
! allocatable arrays
real(8), allocatable :: gvrho(:,:),rfir(:)
complex(8), allocatable :: zfft1(:),zfft2(:)
allocate(gvrho(ngtot,3),rfir(ngtot))
allocate(zfft1(nfgrz),zfft2(nfgrz))
call rzfftifc(3,ngridg,-1,rho,zfft1)
! |∇ρ|
do i=1,3
  do ifg=1,nfgrz
    ig=igrzf(ifg)
    if (ig <= ngvc) then
      zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
    else
      zfft2(ifg)=0.d0
    end if
  end do
  call rzfftifc(3,ngridg,1,gvrho(:,i),zfft2)
end do
grho(:)=sqrt(gvrho(:,1)**2+gvrho(:,2)**2+gvrho(:,3)**2)
! ∇²ρ
do ifg=1,nfgrz
  ig=igrzf(ifg)
  if (ig <= ngvc) then
    zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
  else
    zfft2(ifg)=0.d0
  end if
end do
call rzfftifc(3,ngridg,1,g2rho,zfft2)
! (∇ρ)⋅(∇|∇ρ|)
call rzfftifc(3,ngridg,-1,grho,zfft1)
g3rho(:)=0.d0
do i=1,3
  do ifg=1,nfgrz
    ig=igrzf(ifg)
    zfft2(ifg)=vgc(i,ig)*cmplx(-zfft1(ifg)%im,zfft1(ifg)%re,8)
  end do
  call rzfftifc(3,ngridg,1,rfir,zfft2)
  g3rho(:)=g3rho(:)+gvrho(:,i)*rfir(:)
end do
deallocate(gvrho,rfir,zfft1,zfft2)
end subroutine
!EOC