File: ggamt_2b.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (65 lines) | stat: -rw-r--r-- 1,814 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

! 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: ggamt_2b
! !INTERFACE:
subroutine ggamt_2b(is,g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
! !USES:
use modmain
! !DESCRIPTION:
!   Spin-unpolarised version of {\tt ggamt\_sp\_2b}.
!
! !REVISION HISTORY:
!   Created November 2009 (JKD and TMcQ)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: is
real(8), intent(in) :: g2rho(npmtmax),gvrho(npmtmax,3)
real(8), intent(inout) :: vx(npmtmax),vc(npmtmax)
real(8), intent(in) :: dxdgr2(npmtmax),dcdgr2(npmtmax)
! local variables
integer nr,nri,np,i
! allocatable arrays
real(8), allocatable :: rfmt1(:),rfmt2(:),grfmt(:,:)
allocate(rfmt1(npmtmax),rfmt2(npmtmax),grfmt(npmtmax,3))
nr=nrmt(is)
nri=nrmti(is)
np=npmt(is)
!------------------!
!     exchange     !
!------------------!
! convert dxdgr2 to spherical harmonics
call rfsht(nr,nri,dxdgr2,rfmt1)
! compute grad dxdgr2
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
! (grad dxdgr2).(grad rho) in spherical coordinates
rfmt1(1:np)=0.d0
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),rfmt2)
  rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvrho(1:np,i)
end do
vx(1:np)=vx(1:np)-2.d0*(rfmt1(1:np)+dxdgr2(1:np)*g2rho(1:np))
!---------------------!
!     correlation     !
!---------------------!
! convert dcdgr2 to spherical harmonics
call rfsht(nr,nri,dcdgr2,rfmt1)
! compute grad dcdgr2
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
! (grad dcdgr2).(grad rho) in spherical coordinates
rfmt1(1:np)=0.d0
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),rfmt2)
  rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvrho(1:np,i)
end do
vc(1:np)=vc(1:np)-2.d0*(rfmt1(1:np)+dcdgr2(1:np)*g2rho(1:np))
deallocate(rfmt1,rfmt2,grfmt)
return
end subroutine
!EOC