File: ggair_sp_2a.f90

package info (click to toggle)
elkcode 2.3.22-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 8,972 kB
  • ctags: 3,389
  • sloc: f90: 40,548; fortran: 20,560; perl: 965; makefile: 320; sh: 287; ansic: 67; python: 41
file content (96 lines) | stat: -rw-r--r-- 2,537 bytes parent folder | download | duplicates (4)
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
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96

! Copyright (C) 2009 T. McQueen and J. K. Dewhurst.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: ggair_sp_2a
! !INTERFACE:
subroutine ggair_sp_2a(rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
! !USES:
use modmain
! !DESCRIPTION:
!   Computes the interstitial gradients $\nabla^2\rho^{\uparrow}$,
!   $\nabla^2\rho^{\downarrow}$, $\nabla\rho^{\uparrow}$,
!   $\nabla\rho^{\downarrow}$, $(\nabla\rho^{\uparrow})^2$,
!   $(\nabla\rho^{\downarrow})^2$ and
!   $\nabla\rho^{\uparrow}\cdot\nabla\rho^{\downarrow}$. These are used for GGA
!   functionals of type 2 and meta-GGA. See {\tt ggamt\_sp\_2a} for details.
!
! !REVISION HISTORY:
!   Created November 2009 (JKD and TMcQ)
!EOP
!BOC
implicit none
! arguments
real(8), intent(in) :: rhoup(ngtot)
real(8), intent(in) :: rhodn(ngtot)
real(8), intent(out) :: g2up(ngtot)
real(8), intent(out) :: g2dn(ngtot)
real(8), intent(out) :: gvup(ngtot,3)
real(8), intent(out) :: gvdn(ngtot,3)
real(8), intent(out) :: gup2(ngtot)
real(8), intent(out) :: gdn2(ngtot)
real(8), intent(out) :: gupdn(ngtot)
! local variables
integer ig,ifg,i
! allocatable arrays
complex(8), allocatable :: zfft1(:),zfft2(:)
allocate(zfft1(ngtot),zfft2(ngtot))
!----------------!
!     rho up     !
!----------------!
zfft1(:)=rhoup(:)
call zfftifc(3,ngridg,-1,zfft1)
! compute grad^2 rhoup
zfft2(:)=0.d0
do ig=1,ngvec
  ifg=igfft(ig)
  zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
end do
call zfftifc(3,ngridg,1,zfft2)
g2up(:)=dble(zfft2(:))
! grad rhoup
do i=1,3
  zfft2(:)=0.d0
  do ig=1,ngvec
    ifg=igfft(ig)
    zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8)
  end do
  call zfftifc(3,ngridg,1,zfft2)
  gvup(:,i)=dble(zfft2(:))
end do
! (grad rhoup)^2
gup2(:)=gvup(:,1)**2+gvup(:,2)**2+gvup(:,3)**2
!------------------!
!     rho down     !
!------------------!
zfft1(:)=rhodn(:)
call zfftifc(3,ngridg,-1,zfft1)
! compute grad^2 rhodn
zfft2(:)=0.d0
do ig=1,ngvec
  ifg=igfft(ig)
  zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
end do
call zfftifc(3,ngridg,1,zfft2)
g2dn(:)=dble(zfft2(:))
! grad rhodn
do i=1,3
  zfft2(:)=0.d0
  do ig=1,ngvec
    ifg=igfft(ig)
    zfft2(ifg)=vgc(i,ig)*cmplx(-aimag(zfft1(ifg)),dble(zfft1(ifg)),8)
  end do
  call zfftifc(3,ngridg,1,zfft2)
  gvdn(:,i)=dble(zfft2(:))
end do
! (grad rhodn)^2
gdn2(:)=gvdn(:,1)**2+gvdn(:,2)**2+gvdn(:,3)**2
! (grad rhoup).(grad rhodn)
gupdn(:)=gvup(:,1)*gvdn(:,1)+gvup(:,2)*gvdn(:,2)+gvup(:,3)*gvdn(:,3)
deallocate(zfft1,zfft2)
return
end subroutine
!EOC