The Elk Code
 
Loading...
Searching...
No Matches
ggamt_3.f90
Go to the documentation of this file.
1
2! Copyright (C) 2023 J. K. Dewhurst and S. Sharma.
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_3(is,np,vx,vc,wx,wc,dtdg2r)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: is,np
11real(8), intent(inout) :: vx(np),vc(np)
12real(8), intent(in) :: wx(np),wc(np)
13real(8), intent(in) :: dtdg2r(np)
14! local variables
15integer nr,nri
16! automatic arrays
17real(8) rfmt1(np),rfmt2(np)
18nr=nrmt(is)
19nri=nrmti(is)
20!------------------!
21! exchange !
22!------------------!
23rfmt1(1:np)=wx(1:np)*dtdg2r(1:np)
24call rfsht(nr,nri,rfmt1,rfmt2)
25call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt2,rfmt1)
26call rbsht(nr,nri,rfmt1,rfmt2)
27vx(1:np)=vx(1:np)+rfmt2(1:np)
28!---------------------!
29! correlation !
30!---------------------!
31rfmt1(1:np)=wc(1:np)*dtdg2r(1:np)
32call rfsht(nr,nri,rfmt1,rfmt2)
33call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt2,rfmt1)
34call rbsht(nr,nri,rfmt1,rfmt2)
35vc(1:np)=vc(1:np)+rfmt2(1:np)
36end subroutine
37
subroutine ggamt_3(is, np, vx, vc, wx, wc, dtdg2r)
Definition ggamt_3.f90:7
subroutine grad2rfmt(nr, nri, ri, ri2, wcr, rfmt, g2rfmt)
Definition grad2rfmt.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