The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine ggamt_sp_1(is,np,rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
10! !USES:
11use 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
41implicit none
42! arguments
43integer, intent(in) :: is,np
44real(8), intent(in) :: rhoup(np),rhodn(np)
45real(8), intent(out) :: grho(np),gup(np),gdn(np)
46real(8), intent(out) :: g2up(np),g2dn(np)
47real(8), intent(out) :: g3rho(np),g3up(np),g3dn(np)
48! local variables
49integer nr,nri,i
50! automatic arrays
51real(8) grfmt(np,3),gvup(np,3),gvdn(np,3)
52real(8) rfmt1(np),rfmt2(np)
53nr=nrmt(is)
54nri=nrmti(is)
55!----------------!
56! rho up !
57!----------------!
58! convert rhoup to spherical harmonics
59call rfsht(nr,nri,rhoup,rfmt1)
60! grad rhoup in spherical coordinates
61call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
62do i=1,3
63 call rbsht(nr,nri,grfmt(:,i),gvup(:,i))
64end do
65! |grad rhoup|
66gup(1:np)=sqrt(gvup(1:np,1)**2+gvup(1:np,2)**2+gvup(1:np,3)**2)
67! grad^2 rhoup in spherical coordinates
68call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
69call rbsht(nr,nri,rfmt2,g2up)
70! (grad rhoup).(grad |grad rhoup|)
71call rfsht(nr,nri,gup,rfmt1)
72call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
73g3up(1:np)=0.d0
74do 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)
77end do
78!------------------!
79! rho down !
80!------------------!
81! convert rhodn to spherical harmonics
82call rfsht(nr,nri,rhodn,rfmt1)
83! grad rhodn in spherical coordinates
84call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
85do i=1,3
86 call rbsht(nr,nri,grfmt(:,i),gvdn(:,i))
87end do
88gdn(1:np)=sqrt(gvdn(1:np,1)**2+gvdn(1:np,2)**2+gvdn(1:np,3)**2)
89! grad^2 rhodn in spherical coordinates
90call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
91call rbsht(nr,nri,rfmt2,g2dn)
92! (grad rhodn).(grad |grad rhodn|)
93call rfsht(nr,nri,gdn,rfmt1)
94call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
95g3dn(1:np)=0.d0
96do 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)
99end do
100!-------------!
101! rho !
102!-------------!
103! |grad rho|
104grho(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! (grad rho).(grad |grad rho|)
108call rfsht(nr,nri,grho,rfmt1)
109call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
110g3rho(1:np)=0.d0
111do 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)
114end do
115end subroutine
116!EOC
117
subroutine ggamt_sp_1(is, np, rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn)
subroutine grad2rfmt(nr, nri, ri, ri2, wcr, rfmt, g2rfmt)
Definition grad2rfmt.f90:10
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