The Elk Code
ggamt_sp_1.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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_1
8 ! !INTERFACE:
9 subroutine ggamt_sp_1(is,np,rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! is : species number (in,integer)
14 ! np : number of muffin-tin points (in,integer)
15 ! rhoup : spin-up density in spherical coordinates (in,real(np))
16 ! rhodn : spin-down density (in,real(np))
17 ! grho : |grad rho| (out,real(np))
18 ! gup : |grad rhoup| (out,real(np))
19 ! gdn : |grad rhodn| (out,real(np))
20 ! g2up : grad^2 rhoup (out,real(np))
21 ! g2dn : grad^2 rhodn (out,real(np))
22 ! g3rho : (grad rho).(grad |grad rho|) (out,real(np))
23 ! g3up : (grad rhoup).(grad |grad rhoup|) (out,real(np))
24 ! g3dn : (grad rhodn).(grad |grad rhodn|) (out,real(np))
25 ! !DESCRIPTION:
26 ! Computes $|\nabla\rho|$, $|\nabla\rho^{\uparrow}|$,
27 ! $|\nabla\rho^{\downarrow}|$, $\nabla^2\rho^{\uparrow}$,
28 ! $\nabla^2\rho^{\downarrow}$, $\nabla\rho\cdot(\nabla|\nabla\rho|)$,
29 ! $\nabla\rho^{\uparrow}\cdot(\nabla|\nabla\rho^{\uparrow}|)$ and
30 ! $\nabla\rho^{\downarrow}\cdot(\nabla|\nabla\rho^{\downarrow}|)$
31 ! for a muffin-tin charge density, as required by the generalised gradient
32 ! approximation functionals of type 1 for spin-polarised densities. The input
33 ! densities and output gradients are in terms of spherical coordinates. See
34 ! routines {\tt potxc} and {\tt modxcifc}.
35 !
36 ! !REVISION HISTORY:
37 ! Created April 2004 (JKD)
38 ! Simplified and improved, October 2009 (JKD)
39 !EOP
40 !BOC
41 implicit none
42 ! arguments
43 integer, intent(in) :: is,np
44 real(8), intent(in) :: rhoup(np),rhodn(np)
45 real(8), intent(out) :: grho(np),gup(np),gdn(np)
46 real(8), intent(out) :: g2up(np),g2dn(np)
47 real(8), intent(out) :: g3rho(np),g3up(np),g3dn(np)
48 ! local variables
49 integer nr,nri,i
50 ! automatic arrays
51 real(8) grfmt(np,3),gvup(np,3),gvdn(np,3)
52 real(8) rfmt1(np),rfmt2(np)
53 nr=nrmt(is)
54 nri=nrmti(is)
55 !------------!
56 ! ρ↑ !
57 !------------!
58 ! convert ρ↑ to spherical harmonics
59 call rfsht(nr,nri,rhoup,rfmt1)
60 ! ∇ρ↑ in spherical coordinates
61 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
62 do i=1,3
63  call rbsht(nr,nri,grfmt(:,i),gvup(:,i))
64 end do
65 ! |∇ρ↑|
66 gup(1:np)=sqrt(gvup(1:np,1)**2+gvup(1:np,2)**2+gvup(1:np,3)**2)
67 ! ∇²ρ↑ in spherical coordinates
68 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
69 call rbsht(nr,nri,rfmt2,g2up)
70 ! (∇ρ↑)⋅(∇|∇ρ↑|)
71 call rfsht(nr,nri,gup,rfmt1)
72 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
73 g3up(1:np)=0.d0
74 do i=1,3
75  call rbsht(nr,nri,grfmt(:,i),rfmt1)
76  g3up(1:np)=g3up(1:np)+gvup(1:np,i)*rfmt1(1:np)
77 end do
78 !------------!
79 ! ρ↓ !
80 !------------!
81 ! convert ρ↓ to spherical harmonics
82 call rfsht(nr,nri,rhodn,rfmt1)
83 ! ∇ρ↓ in spherical coordinates
84 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
85 do i=1,3
86  call rbsht(nr,nri,grfmt(:,i),gvdn(:,i))
87 end do
88 gdn(1:np)=sqrt(gvdn(1:np,1)**2+gvdn(1:np,2)**2+gvdn(1:np,3)**2)
89 ! ∇²ρ↓ in spherical coordinates
90 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
91 call rbsht(nr,nri,rfmt2,g2dn)
92 ! (∇ρ↓)⋅(∇|∇ρ↓|)
93 call rfsht(nr,nri,gdn,rfmt1)
94 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
95 g3dn(1:np)=0.d0
96 do i=1,3
97  call rbsht(nr,nri,grfmt(:,i),rfmt1)
98  g3dn(1:np)=g3dn(1:np)+gvdn(1:np,i)*rfmt1(1:np)
99 end do
100 !-----------!
101 ! ρ !
102 !-----------!
103 ! |∇ρ|
104 grho(1:np)=sqrt((gvup(1:np,1)+gvdn(1:np,1))**2 &
105  +(gvup(1:np,2)+gvdn(1:np,2))**2 &
106  +(gvup(1:np,3)+gvdn(1:np,3))**2)
107 ! (∇ρ)⋅(∇|∇ρ|)
108 call rfsht(nr,nri,grho,rfmt1)
109 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
110 g3rho(1:np)=0.d0
111 do i=1,3
112  call rbsht(nr,nri,grfmt(:,i),rfmt1)
113  g3rho(1:np)=g3rho(1:np)+(gvup(1:np,i)+gvdn(1:np,i))*rfmt1(1:np)
114 end do
115 end subroutine
116 !EOC
117 
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
subroutine ggamt_sp_1(is, np, rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn)
Definition: ggamt_sp_1.f90:10
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
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