The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine ggamt_2b(is,np,g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
10! !USES:
11use 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
19implicit none
20! arguments
21integer, intent(in) :: is,np
22real(8), intent(in) :: g2rho(np),gvrho(np,3)
23real(8), intent(inout) :: vx(np),vc(np)
24real(8), intent(in) :: dxdgr2(np),dcdgr2(np)
25! local variables
26integer nr,nri,i
27! automatic arrays
28real(8) rfmt1(np),rfmt2(np),grfmt(np,3)
29nr=nrmt(is)
30nri=nrmti(is)
31!------------------!
32! exchange !
33!------------------!
34! convert dxdgr2 to spherical harmonics
35call rfsht(nr,nri,dxdgr2,rfmt1)
36! compute grad dxdgr2
37call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
38! (grad dxdgr2).(grad rho) in spherical coordinates
39rfmt1(1:np)=0.d0
40do 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)
43end do
44vx(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
49call rfsht(nr,nri,dcdgr2,rfmt1)
50! compute grad dcdgr2
51call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
52! (grad dcdgr2).(grad rho) in spherical coordinates
53rfmt1(1:np)=0.d0
54do 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)
57end do
58vc(1:np)=vc(1:np)-2.d0*(rfmt1(1:np)+dcdgr2(1:np)*g2rho(1:np))
59end subroutine
60!EOC
61
subroutine ggamt_2b(is, np, g2rho, gvrho, vx, vc, dxdgr2, dcdgr2)
Definition ggamt_2b.f90:10
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition gradrfmt.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition rfsht.f90:7