The Elk Code
ggamt_4.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine ggamt_4(is,np,gvrho,vx,vc,wx,wc,dtdr,dtdgr2)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,np
11 real(8), intent(in) :: gvrho(np,3)
12 real(8), intent(inout) :: vx(np),vc(np)
13 real(8), intent(in) :: wx(np),wc(np)
14 real(8), intent(in) :: dtdr(np),dtdgr2(np)
15 ! local variables
16 integer nr,nri,i
17 ! automatic arrays
18 real(8) grfmt(np,3),rfmt1(np),rfmt2(np),rfmt3(np)
19 nr=nrmt(is)
20 nri=nrmti(is)
21 !------------------!
22 ! exchange !
23 !------------------!
24 vx(1:np)=vx(1:np)+wx(1:np)*dtdr(1:np)
25 rfmt1(1:np)=wx(1:np)*dtdgr2(1:np)
26 do i=1,3
27  rfmt2(1:np)=rfmt1(1:np)*gvrho(1:np,i)
28  call rfsht(nr,nri,rfmt2,rfmt3)
29  call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt3,np,grfmt)
30  call rbsht(nr,nri,grfmt(:,i),rfmt2)
31  vx(1:np)=vx(1:np)-2.d0*rfmt2(1:np)
32 end do
33 !---------------------!
34 ! correlation !
35 !---------------------!
36 vc(1:np)=vc(1:np)+wc(1:np)*dtdr(1:np)
37 rfmt1(1:np)=wc(1:np)*dtdgr2(1:np)
38 do i=1,3
39  rfmt2(1:np)=rfmt1(1:np)*gvrho(1:np,i)
40  call rfsht(nr,nri,rfmt2,rfmt3)
41  call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt3,np,grfmt)
42  call rbsht(nr,nri,grfmt(:,i),rfmt2)
43  vc(1:np)=vc(1:np)-2.d0*rfmt2(1:np)
44 end do
45 end subroutine
46 
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
subroutine ggamt_4(is, np, gvrho, vx, vc, wx, wc, dtdr, dtdgr2)
Definition: ggamt_4.f90:7
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition: gradrfmt.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