The Elk Code
 
Loading...
Searching...
No Matches
ggamt_sp_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_sp_2b
8! !INTERFACE:
9subroutine ggamt_sp_2b(is,np,g2up,g2dn,gvup,gvdn,vxup,vxdn,vcup,vcdn,dxdgu2, &
10 dxdgd2,dxdgud,dcdgu2,dcdgd2,dcdgud)
11! !USES:
12use modmain
13! !DESCRIPTION:
14! Post processing step of muffin-tin gradients for GGA type 2. See routine
15! {\tt ggamt\_sp\_2a} for full details.
16!
17! !REVISION HISTORY:
18! Created November 2009 (JKD and TMcQ)
19!EOP
20!BOC
21implicit none
22! arguments
23integer, intent(in) :: is,np
24real(8), intent(in) :: g2up(np),g2dn(np)
25real(8), intent(in) :: gvup(np,3),gvdn(np,3)
26real(8), intent(inout) :: vxup(np),vxdn(np)
27real(8), intent(inout) :: vcup(np),vcdn(np)
28real(8), intent(in) :: dxdgu2(np),dxdgd2(np),dxdgud(np)
29real(8), intent(in) :: dcdgu2(np),dcdgd2(np),dcdgud(np)
30! local variables
31integer nr,nri,i
32! automatic arrays
33real(8) rfmt1(np),rfmt2(np),grfmt(np,3)
34nr=nrmt(is)
35nri=nrmti(is)
36!------------------!
37! exchange !
38!------------------!
39! convert dxdgu2 to spherical harmonics
40call rfsht(nr,nri,dxdgu2,rfmt1)
41! compute grad dxdgu2
42call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
43! (grad dxdgu2).(grad rhoup) in spherical coordinates
44rfmt1(1:np)=0.d0
45do i=1,3
46 call rbsht(nr,nri,grfmt(:,i),rfmt2)
47 rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvup(1:np,i)
48end do
49vxup(1:np)=vxup(1:np)-2.d0*(rfmt1(1:np)+dxdgu2(1:np)*g2up(1:np)) &
50 -dxdgud(1:np)*g2dn(1:np)
51! convert dxdgd2 to spherical harmonics
52call rfsht(nr,nri,dxdgd2,rfmt1)
53! compute grad dxdgd2
54call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
55! (grad dxdgd2).(grad rhodn) in spherical coordinates
56rfmt1(1:np)=0.d0
57do i=1,3
58 call rbsht(nr,nri,grfmt(:,i),rfmt2)
59 rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvdn(1:np,i)
60end do
61vxdn(1:np)=vxdn(1:np)-2.d0*(rfmt1(1:np)+dxdgd2(1:np)*g2dn(1:np)) &
62 -dxdgud(1:np)*g2up(1:np)
63! convert dxdgud to spherical harmonics
64call rfsht(nr,nri,dxdgud,rfmt1)
65! compute grad dxdgud
66call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
67! (grad dxdgud).(grad rhodn) and (grad dxdgud).(grad rhoup)
68do i=1,3
69 call rbsht(nr,nri,grfmt(:,i),rfmt1)
70 vxup(1:np)=vxup(1:np)-rfmt1(1:np)*gvdn(1:np,i)
71 vxdn(1:np)=vxdn(1:np)-rfmt1(1:np)*gvup(1:np,i)
72end do
73!---------------------!
74! correlation !
75!---------------------!
76! convert dcdgu2 to spherical harmonics
77call rfsht(nr,nri,dcdgu2,rfmt1)
78! compute grad dcdgu2
79call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
80! (grad dcdgu2).(grad rhoup) in spherical coordinates
81rfmt1(1:np)=0.d0
82do i=1,3
83 call rbsht(nr,nri,grfmt(:,i),rfmt2)
84 rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvup(1:np,i)
85end do
86vcup(1:np)=vcup(1:np)-2.d0*(rfmt1(1:np)+dcdgu2(1:np)*g2up(1:np)) &
87 -dcdgud(1:np)*g2dn(1:np)
88! convert dcdgd2 to spherical harmonics
89call rfsht(nr,nri,dcdgd2,rfmt1)
90! compute grad dcdgd2
91call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
92! (grad dcdgd2).(grad rhodn) in spherical coordinates
93rfmt1(1:np)=0.d0
94do i=1,3
95 call rbsht(nr,nri,grfmt(:,i),rfmt2)
96 rfmt1(1:np)=rfmt1(1:np)+rfmt2(1:np)*gvdn(1:np,i)
97end do
98vcdn(1:np)=vcdn(1:np)-2.d0*(rfmt1(1:np)+dcdgd2(1:np)*g2dn(1:np)) &
99 -dcdgud(1:np)*g2up(1:np)
100! convert dcdgud to spherical harmonics
101call rfsht(nr,nri,dcdgud,rfmt1)
102! compute grad dcdgud
103call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
104! (grad dcdgud).(grad rhodn) and (grad dcdgud).(grad rhoup)
105do i=1,3
106 call rbsht(nr,nri,grfmt(:,i),rfmt1)
107 vcup(1:np)=vcup(1:np)-rfmt1(1:np)*gvdn(1:np,i)
108 vcdn(1:np)=vcdn(1:np)-rfmt1(1:np)*gvup(1:np,i)
109end do
110end subroutine
111!EOC
112
subroutine ggamt_sp_2b(is, np, g2up, g2dn, gvup, gvdn, vxup, vxdn, vcup, vcdn, dxdgu2, dxdgd2, dxdgud, dcdgu2, dcdgd2, dcdgud)
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