The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine ggamt_4(is,np,gvrho,vx,vc,wx,wc,dtdr,dtdgr2)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: is,np
11real(8), intent(in) :: gvrho(np,3)
12real(8), intent(inout) :: vx(np),vc(np)
13real(8), intent(in) :: wx(np),wc(np)
14real(8), intent(in) :: dtdr(np),dtdgr2(np)
15! local variables
16integer nr,nri,i
17! automatic arrays
18real(8) grfmt(np,3),rfmt1(np),rfmt2(np),rfmt3(np)
19nr=nrmt(is)
20nri=nrmti(is)
21!------------------!
22! exchange !
23!------------------!
24vx(1:np)=vx(1:np)+wx(1:np)*dtdr(1:np)
25rfmt1(1:np)=wx(1:np)*dtdgr2(1:np)
26do 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)
32end do
33!---------------------!
34! correlation !
35!---------------------!
36vc(1:np)=vc(1:np)+wc(1:np)*dtdr(1:np)
37rfmt1(1:np)=wc(1:np)*dtdgr2(1:np)
38do 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)
44end do
45end subroutine
46
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
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