The Elk Code
ggamt_sp_2a.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_2a
8 ! !INTERFACE:
9 subroutine ggamt_sp_2a(is,np,rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Computes the muffin-tin gradients $\nabla^2\rho^{\uparrow}$,
14 ! $\nabla^2\rho^{\downarrow}$, $\nabla\rho^{\uparrow}$,
15 ! $\nabla\rho^{\downarrow}$, $(\nabla\rho^{\uparrow})^2$,
16 ! $(\nabla\rho^{\downarrow})^2$ and
17 ! $\nabla\rho^{\uparrow}\cdot\nabla\rho^{\downarrow}$, which are passed in to
18 ! GGA functional subroutines of type 2. The exchange-correlation energy in
19 ! these routines has the functional form
20 ! $$ E_{xc}[\rho^{\uparrow},\rho^{\downarrow}]=\int d^3r\,\hat{\epsilon}_{xc}
21 ! \bigl(\rho^{\uparrow}({\bf r}),\rho^{\downarrow}({\bf r}),
22 ! (\nabla\rho^{\uparrow}({\bf r}))^2,(\nabla\rho^{\downarrow}({\bf r}))^2,
23 ! \nabla\rho^{\uparrow}({\bf r})
24 ! \cdot\nabla\rho^{\downarrow}({\bf r})\bigr), $$
25 ! where $\hat{\epsilon}_{xc}({\bf r})=\epsilon_{xc}({\bf r})\rho({\bf r})$ is
26 ! the xc energy per unit volume, with $\epsilon_{xc}$ being the xc energy per
27 ! electron, and $\rho=\rho^{\uparrow}+\rho^{\downarrow}$. From the gradients
28 ! above, type 2 GGA routines return $\epsilon_{xc}$, but not directly the xc
29 ! potentials. Instead they generate the derivatives
30 ! $\partial\hat{\epsilon}_{xc}/\partial\rho^{\uparrow}({\bf r})$,
31 ! $\partial\hat{\epsilon}_{xc}/\partial(\nabla\rho^{\uparrow}({\bf r}))^2$,
32 ! and the same for down spin, as well as
33 ! $\partial\hat{\epsilon}_{xc}/\partial(\nabla\rho^{\uparrow}({\bf r})
34 ! \cdot\nabla\rho^{\downarrow}({\bf r}))$. In a post-processing step invoked
35 ! by {\tt ggamt\_sp\_2b}, integration by parts is used to obtain the xc
36 ! potential explicitly with
37 ! \begin{align*}
38 ! V_{xc}^{\uparrow}({\bf r})=&\frac{\partial\hat{\epsilon}_{xc}}{\partial
39 ! \rho^{\uparrow}({\bf r})}-2\left(\nabla\frac{\partial\hat{\epsilon}_{xc}}
40 ! {\partial(\nabla\rho^{\uparrow})^2}\right)\cdot\nabla\rho^{\uparrow}
41 ! -2\frac{\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow})^2}\nabla^2
42 ! \rho^{\uparrow}\\
43 ! &-\left(\nabla\frac{\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow}
44 ! \cdot\nabla\rho^{\downarrow})}\right)\cdot\nabla\rho^{\downarrow}
45 ! -\frac{\partial\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow}\cdot
46 ! \nabla\rho^{\downarrow})}\nabla^2\rho^{\downarrow},
47 ! \end{align*}
48 ! and similarly for $V_{xc}^{\downarrow}$.
49 !
50 ! !REVISION HISTORY:
51 ! Created November 2009 (JKD and TMcQ)
52 !EOP
53 !BOC
54 implicit none
55 ! arguments
56 integer, intent(in) :: is,np
57 real(8), intent(in) :: rhoup(np),rhodn(np)
58 real(8), intent(out) :: g2up(np),g2dn(np)
59 real(8), intent(out) :: gvup(np,3),gvdn(np,3)
60 real(8), intent(out) :: gup2(np),gdn2(np),gupdn(np)
61 ! local variables
62 integer nr,nri,i
63 ! automatic arrays
64 real(8) rfmt1(np),rfmt2(np),grfmt(np,3)
65 nr=nrmt(is)
66 nri=nrmti(is)
67 !------------!
68 ! ρ↑ !
69 !------------!
70 ! convert ρ↑ to spherical harmonics
71 call rfsht(nr,nri,rhoup,rfmt1)
72 ! compute ∇²ρ↑ in spherical coordinates
73 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
74 call rbsht(nr,nri,rfmt2,g2up)
75 ! ∇ρ↑ in spherical coordinates
76 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
77 do i=1,3
78  call rbsht(nr,nri,grfmt(:,i),gvup(:,i))
79 end do
80 ! (∇ρ↑)²
81 gup2(1:np)=gvup(1:np,1)**2+gvup(1:np,2)**2+gvup(1:np,3)**2
82 !------------!
83 ! ρ↓ !
84 !------------!
85 ! convert ρ↓ to spherical harmonics
86 call rfsht(nr,nri,rhodn,rfmt1)
87 ! compute ∇²ρ↓ in spherical coordinates
88 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
89 call rbsht(nr,nri,rfmt2,g2dn)
90 ! ∇ρ↓ in spherical coordinates
91 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
92 do i=1,3
93  call rbsht(nr,nri,grfmt(:,i),gvdn(:,i))
94 end do
95 ! (∇ρ↓)²
96 gdn2(1:np)=gvdn(1:np,1)**2+gvdn(1:np,2)**2+gvdn(1:np,3)**2
97 ! (∇ρ↑)⋅(∇ρ↓)
98 gupdn(1:np)=gvup(1:np,1)*gvdn(1:np,1) &
99  +gvup(1:np,2)*gvdn(1:np,2) &
100  +gvup(1:np,3)*gvdn(1:np,3)
101 end subroutine
102 !EOC
103 
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.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
subroutine grad2rfmt(nr, nri, ri, ri2, wcr, rfmt, g2rfmt)
Definition: grad2rfmt.f90:10
subroutine ggamt_sp_2a(is, np, rhoup, rhodn, g2up, g2dn, gvup, gvdn, gup2, gdn2, gupdn)
Definition: ggamt_sp_2a.f90:10
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