File: ggamt_sp_1.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 (121 lines) | stat: -rw-r--r-- 4,177 bytes parent folder | download | duplicates (2)
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
108
109
110
111
112
113
114
115
116
117
118
119
120
121

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: ggamt_sp_1
! !INTERFACE:
subroutine ggamt_sp_1(is,rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
! !USES:
use modmain
! !INPUT/OUTPUT PARAMETERS:
!   is    : species number (in,integer)
!   rhoup : spin-up density in spherical coordinates (in,real(npmtmax))
!   rhodn : spin-down density (in,real(npmtmax))
!   grho  : |grad rho| (out,real(npmtmax))
!   gup   : |grad rhoup| (out,real(npmtmax))
!   gdn   : |grad rhodn| (out,real(npmtmax))
!   g2up  : grad^2 rhoup (out,real(npmtmax))
!   g2dn  : grad^2 rhodn (out,real(npmtmax))
!   g3rho : (grad rho).(grad |grad rho|) (out,real(npmtmax))
!   g3up  : (grad rhoup).(grad |grad rhoup|) (out,real(npmtmax))
!   g3dn  : (grad rhodn).(grad |grad rhodn|) (out,real(npmtmax))
! !DESCRIPTION:
!   Computes $|\nabla\rho|$, $|\nabla\rho^{\uparrow}|$,
!   $|\nabla\rho^{\downarrow}|$, $\nabla^2\rho^{\uparrow}$,
!   $\nabla^2\rho^{\downarrow}$, $\nabla\rho\cdot(\nabla|\nabla\rho|)$,
!   $\nabla\rho^{\uparrow}\cdot(\nabla|\nabla\rho^{\uparrow}|)$ and
!   $\nabla\rho^{\downarrow}\cdot(\nabla|\nabla\rho^{\downarrow}|)$
!   for a muffin-tin charge density, as required by the generalised gradient
!   approximation functionals of type 1 for spin-polarised densities. The input
!   densities and output gradients are in terms of spherical coordinates. See
!   routines {\tt potxc} and {\tt modxcifc}.
!
! !REVISION HISTORY:
!   Created April 2004 (JKD)
!   Simplified and improved, October 2009 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: is
real(8), intent(in) :: rhoup(npmtmax),rhodn(npmtmax)
real(8), intent(out) :: grho(npmtmax),gup(npmtmax),gdn(npmtmax)
real(8), intent(out) :: g2up(npmtmax),g2dn(npmtmax)
real(8), intent(out) :: g3rho(npmtmax),g3up(npmtmax),g3dn(npmtmax)
! local variables
integer nr,nri,np,i
! allocatable arrays
real(8), allocatable :: rfmt1(:),rfmt2(:),grfmt(:,:)
real(8), allocatable :: gvup(:,:),gvdn(:,:)
allocate(rfmt1(npmtmax),rfmt2(npmtmax))
allocate(grfmt(npmtmax,3),gvup(npmtmax,3),gvdn(npmtmax,3))
nr=nrmt(is)
nri=nrmti(is)
np=npmt(is)
!----------------!
!     rho up     !
!----------------!
! convert rhoup to spherical harmonics
call rfsht(nr,nri,rhoup,rfmt1)
! 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|
gup(1:np)=sqrt(gvup(1:np,1)**2+gvup(1:np,2)**2+gvup(1:np,3)**2)
! grad^2 rhoup in spherical coordinates
call grad2rfmt(nr,nri,rsp(:,is),rfmt1,rfmt2)
call rbsht(nr,nri,rfmt2,g2up)
! (grad rhoup).(grad |grad rhoup|)
call rfsht(nr,nri,gup,rfmt1)
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
g3up(1:np)=0.d0
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),rfmt1)
  g3up(1:np)=g3up(1:np)+gvup(1:np,i)*rfmt1(1:np)
end do
!------------------!
!     rho down     !
!------------------!
! convert rhodn to spherical harmonics
call rfsht(nr,nri,rhodn,rfmt1)
! 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
gdn(1:np)=sqrt(gvdn(1:np,1)**2+gvdn(1:np,2)**2+gvdn(1:np,3)**2)
! grad^2 rhodn in spherical coordinates
call grad2rfmt(nr,nri,rsp(:,is),rfmt1,rfmt2)
call rbsht(nr,nri,rfmt2,g2dn)
! (grad rhodn).(grad |grad rhodn|)
call rfsht(nr,nri,gdn,rfmt1)
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
g3dn(1:np)=0.d0
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),rfmt1)
  g3dn(1:np)=g3dn(1:np)+gvdn(1:np,i)*rfmt1(1:np)
end do
!-------------!
!     rho     !
!-------------!
! |grad rho|
grho(1:np)=sqrt((gvup(1:np,1)+gvdn(1:np,1))**2 &
               +(gvup(1:np,2)+gvdn(1:np,2))**2 &
               +(gvup(1:np,3)+gvdn(1:np,3))**2)
! (grad rho).(grad |grad rho|)
call rfsht(nr,nri,grho,rfmt1)
call gradrfmt(nr,nri,rsp(:,is),rfmt1,npmtmax,grfmt)
g3rho(1:np)=0.d0
do i=1,3
  call rbsht(nr,nri,grfmt(:,i),rfmt1)
  g3rho(1:np)=g3rho(1:np)+(gvup(1:np,i)+gvdn(1:np,i))*rfmt1(1:np)
end do
deallocate(rfmt1,rfmt2,grfmt,gvup,gvdn)
return
end subroutine
!EOC