File: ggamt_sp_2a.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (107 lines) | stat: -rw-r--r-- 4,292 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107

! Copyright (C) 2009 T. McQueen and J. K. Dewhurst.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: ggamt_sp_2a
! !INTERFACE:
subroutine ggamt_sp_2a(is,rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
! !USES:
use modmain
! !DESCRIPTION:
!   Computes the muffin-tin gradients $\nabla^2\rho^{\uparrow}$,
!   $\nabla^2\rho^{\downarrow}$, $\nabla\rho^{\uparrow}$,
!   $\nabla\rho^{\downarrow}$, $(\nabla\rho^{\uparrow})^2$,
!   $(\nabla\rho^{\downarrow})^2$ and
!   $\nabla\rho^{\uparrow}\cdot\nabla\rho^{\downarrow}$, which are passed in to
!   GGA functional subroutines of type 2. The exchange-correlation energy in
!   these routines has the functional form
!   $$ E_{xc}[\rho^{\uparrow},\rho^{\downarrow}]=\int d^3r\,\hat{\epsilon}_{xc}
!    \bigl(\rho^{\uparrow}({\bf r}),\rho^{\downarrow}({\bf r}),
!    (\nabla\rho^{\uparrow}({\bf r}))^2,(\nabla\rho^{\downarrow}({\bf r}))^2,
!    \nabla\rho^{\uparrow}({\bf r})
!    \cdot\nabla\rho^{\downarrow}({\bf r})\bigr), $$
!   where $\hat{\epsilon}_{xc}({\bf r})=\epsilon_{xc}({\bf r})\rho({\bf r})$ is
!   the xc energy per unit volume, with $\epsilon_{xc}$ being the xc energy per
!   electron, and $\rho=\rho^{\uparrow}+\rho^{\downarrow}$. From the gradients
!   above, type 2 GGA routines return $\epsilon_{xc}$, but not directly the xc
!   potentials. Instead they generate the derivatives
!   $\partial\hat{\epsilon}_{xc}/\partial\rho^{\uparrow}({\bf r})$,
!   $\partial\hat{\epsilon}_{xc}/\partial(\nabla\rho^{\uparrow}({\bf r}))^2$,
!   and the same for down spin, as well as
!   $\partial\hat{\epsilon}_{xc}/\partial(\nabla\rho^{\uparrow}({\bf r})
!   \cdot\nabla\rho^{\downarrow}({\bf r}))$. In a post-processing step invoked
!   by {\tt ggamt\_sp\_2b}, integration by parts is used to obtain the xc
!   potential explicitly with
!   \begin{align*}
!    V_{xc}^{\uparrow}({\bf r})=&\frac{\partial\hat{\epsilon}_{xc}}{\partial
!    \rho^{\uparrow}({\bf r})}-2\left(\nabla\frac{\partial\hat{\epsilon}_{xc}}
!    {\partial(\nabla\rho^{\uparrow})^2}\right)\cdot\nabla\rho^{\uparrow}
!    -2\frac{\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow})^2}\nabla^2
!    \rho^{\uparrow}\\
!    &-\left(\nabla\frac{\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow}
!    \cdot\nabla\rho^{\downarrow})}\right)\cdot\nabla\rho^{\downarrow}
!    -\frac{\partial\hat{\epsilon}_{xc}}{\partial(\nabla\rho^{\uparrow}\cdot
!    \nabla\rho^{\downarrow})}\nabla^2\rho^{\downarrow},
!   \end{align*}
!   and similarly for $V_{xc}^{\downarrow}$.
!
! !REVISION HISTORY:
!   Created November 2009 (JKD and TMcQ)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: is
real(8), intent(in) :: rhoup(npmtmax),rhodn(npmtmax)
real(8), intent(out) :: g2up(npmtmax),g2dn(npmtmax)
real(8), intent(out) :: gvup(npmtmax,3),gvdn(npmtmax,3)
real(8), intent(out) :: gup2(npmtmax),gdn2(npmtmax),gupdn(npmtmax)
! local variables
integer nr,nri,np,i
! allocatable arrays
real(8), allocatable :: rfmt1(:),rfmt2(:),grfmt(:,:)
allocate(rfmt1(npmtmax),rfmt2(npmtmax),grfmt(npmtmax,3))
nr=nrmt(is)
nri=nrmti(is)
np=npmt(is)
!----------------!
!     rho up     !
!----------------!
! convert rhoup to spherical harmonics
call rfsht(nr,nri,rhoup,rfmt1)
! compute grad^2 rhoup in spherical coordinates
call grad2rfmt(nr,nri,rsp(:,is),rfmt1,rfmt2)
call rbsht(nr,nri,rfmt2,g2up)
! grad rhoup in spherical coordinates
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),gvup(:,i))
end do
! (grad rhoup)^2
gup2(1:np)=gvup(1:np,1)**2+gvup(1:np,2)**2+gvup(1:np,3)**2
!------------------!
!     rho down     !
!------------------!
! convert rhodn to spherical harmonics
call rfsht(nr,nri,rhodn,rfmt1)
! compute grad^2 rhodn in spherical coordinates
call grad2rfmt(nr,nri,rsp(:,is),rfmt1,rfmt2)
call rbsht(nr,nri,rfmt2,g2dn)
! grad rhodn in spherical coordinates
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),gvdn(:,i))
end do
! (grad rhodn)^2
gdn2(1:np)=gvdn(1:np,1)**2+gvdn(1:np,2)**2+gvdn(1:np,3)**2
! (grad rhoup).(grad rhodn)
gupdn(1:np)=gvup(1:np,1)*gvdn(1:np,1) &
           +gvup(1:np,2)*gvdn(1:np,2) &
           +gvup(1:np,3)*gvdn(1:np,3)
deallocate(rfmt1,rfmt2,grfmt)
return
end subroutine
!EOC