The Elk Code
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:
9 subroutine ggamt_sp_2b(is,np,g2up,g2dn,gvup,gvdn,vxup,vxdn,vcup,vcdn,dxdgu2, &
10  dxdgd2,dxdgud,dcdgu2,dcdgd2,dcdgud)
11 ! !USES:
12 use 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
21 implicit none
22 ! arguments
23 integer, intent(in) :: is,np
24 real(8), intent(in) :: g2up(np),g2dn(np)
25 real(8), intent(in) :: gvup(np,3),gvdn(np,3)
26 real(8), intent(inout) :: vxup(np),vxdn(np)
27 real(8), intent(inout) :: vcup(np),vcdn(np)
28 real(8), intent(in) :: dxdgu2(np),dxdgd2(np),dxdgud(np)
29 real(8), intent(in) :: dcdgu2(np),dcdgd2(np),dcdgud(np)
30 ! local variables
31 integer nr,nri,i
32 ! automatic arrays
33 real(8) rfmt1(np),rfmt2(np),grfmt(np,3)
34 nr=nrmt(is)
35 nri=nrmti(is)
36 !------------------!
37 ! exchange !
38 !------------------!
39 ! convert dxdgu2 to spherical harmonics
40 call rfsht(nr,nri,dxdgu2,rfmt1)
41 ! compute ∇ dxdgu2
42 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
43 ! (∇ dxdgu2)⋅(∇ρ↑) in spherical coordinates
44 rfmt1(1:np)=0.d0
45 do 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)
48 end do
49 vxup(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
52 call rfsht(nr,nri,dxdgd2,rfmt1)
53 ! compute ∇ dxdgd2
54 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
55 ! (∇ dxdgd2)⋅(∇ρ↓) in spherical coordinates
56 rfmt1(1:np)=0.d0
57 do 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)
60 end do
61 vxdn(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
64 call rfsht(nr,nri,dxdgud,rfmt1)
65 ! compute ∇ dxdgud
66 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
67 ! (∇ dxdgud)⋅(∇ρ↓) and (∇ dxdgud)⋅(∇ρ↑)
68 do 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)
72 end do
73 !---------------------!
74 ! correlation !
75 !---------------------!
76 ! convert dcdgu2 to spherical harmonics
77 call rfsht(nr,nri,dcdgu2,rfmt1)
78 ! compute ∇ dcdgu2
79 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
80 ! (∇ dcdgu2)⋅(∇ρ↑) in spherical coordinates
81 rfmt1(1:np)=0.d0
82 do 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)
85 end do
86 vcup(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
89 call rfsht(nr,nri,dcdgd2,rfmt1)
90 ! compute ∇ dcdgd2
91 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
92 ! (∇ dcdgd2)⋅(∇ρ↓) in spherical coordinates
93 rfmt1(1:np)=0.d0
94 do 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)
97 end do
98 vcdn(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
101 call rfsht(nr,nri,dcdgud,rfmt1)
102 ! compute ∇ dcdgud
103 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
104 ! (∇ dcdgud)⋅(∇ρ↓) and (∇ dcdgud)⋅(∇ρ↑)
105 do 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)
109 end do
110 end subroutine
111 !EOC
112 
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
subroutine ggamt_sp_2b(is, np, g2up, g2dn, gvup, gvdn, vxup, vxdn, vcup, vcdn, dxdgu2, dxdgd2, dxdgud, dcdgu2, dcdgd2, dcdgud)
Definition: ggamt_sp_2b.f90:11
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