The Elk Code
ggamt_2b.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: ggamt_2b
8 ! !INTERFACE:
9 subroutine ggamt_2b(is,np,g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Spin-unpolarised version of {\tt ggamt\_sp\_2b}.
14 !
15 ! !REVISION HISTORY:
16 ! Created November 2009 (JKD and TMcQ)
17 !EOP
18 !BOC
19 implicit none
20 ! arguments
21 integer, intent(in) :: is,np
22 real(8), intent(in) :: g2rho(np),gvrho(np,3)
23 real(8), intent(inout) :: vx(np),vc(np)
24 real(8), intent(in) :: dxdgr2(np),dcdgr2(np)
25 ! local variables
26 integer nr,nri,i
27 ! automatic arrays
28 real(8) rfmt1(np),rfmt2(np),grfmt(np,3)
29 nr=nrmt(is)
30 nri=nrmti(is)
31 !------------------!
32 ! exchange !
33 !------------------!
34 ! convert dxdgr2 to spherical harmonics
35 call rfsht(nr,nri,dxdgr2,rfmt1)
36 ! compute ∇ dxdgr2
37 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
38 ! (∇ dxdgr2)⋅(∇ρ) in spherical coordinates
39 rfmt1(1:np)=0.d0
40 do i=1,3
41  call rbsht(nr,nri,grfmt(:,i),rfmt2)
42  rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvrho(1:np,i)
43 end do
44 vx(1:np)=vx(1:np)-2.d0*(rfmt1(1:np)+dxdgr2(1:np)*g2rho(1:np))
45 !---------------------!
46 ! correlation !
47 !---------------------!
48 ! convert dcdgr2 to spherical harmonics
49 call rfsht(nr,nri,dcdgr2,rfmt1)
50 ! compute ∇ dcdgr2
51 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
52 ! (∇ dcdgr2)⋅(∇ρ) in spherical coordinates
53 rfmt1(1:np)=0.d0
54 do i=1,3
55  call rbsht(nr,nri,grfmt(:,i),rfmt2)
56  rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvrho(1:np,i)
57 end do
58 vc(1:np)=vc(1:np)-2.d0*(rfmt1(1:np)+dcdgr2(1:np)*g2rho(1:np))
59 end subroutine
60 !EOC
61 
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition: gradrfmt.f90:10
subroutine ggamt_2b(is, np, g2rho, gvrho, vx, vc, dxdgr2, dcdgr2)
Definition: ggamt_2b.f90:10
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition: rfsht.f90:7
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150