The Elk Code
 
Loading...
Searching...
No Matches
ggamt_1.f90
Go to the documentation of this file.
1
2! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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_1
8! !INTERFACE:
9subroutine ggamt_1(tsh,is,np,rho,grho,g2rho,g3rho)
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Spin-unpolarised version of {\tt ggamt\_sp\_1}.
14!
15! !REVISION HISTORY:
16! Created November 2009 (JKD)
17!EOP
18!BOC
19implicit none
20! arguments
21logical, intent(in) :: tsh
22integer, intent(in) :: is,np
23real(8), intent(in) :: rho(np)
24real(8), intent(out) :: grho(np),g2rho(np),g3rho(np)
25! local variables
26integer nr,nri,i
27! automatic arrays
28real(8) grfmt(np,3),gvrho(np,3),rfmt1(np),rfmt2(np)
29nr=nrmt(is)
30nri=nrmti(is)
31! |grad rho|
32if (tsh) then
33 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rho,np,grfmt)
34else
35 call rfsht(nr,nri,rho,rfmt1)
36 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
37end if
38do i=1,3
39 call rbsht(nr,nri,grfmt(:,i),gvrho(:,i))
40end do
41grho(1:np)=sqrt(gvrho(1:np,1)**2+gvrho(1:np,2)**2+gvrho(1:np,3)**2)
42! grad^2 rho in spherical coordinates
43if (tsh) then
44 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rho,rfmt2)
45else
46 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
47end if
48call rbsht(nr,nri,rfmt2,g2rho)
49! (grad rho).(grad |grad rho|)
50call rfsht(nr,nri,grho,rfmt2)
51call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt2,np,grfmt)
52g3rho(1:np)=0.d0
53do i=1,3
54 call rbsht(nr,nri,grfmt(:,i),rfmt2)
55 g3rho(1:np)=g3rho(1:np)+gvrho(1:np,i)*rfmt2(1:np)
56end do
57end subroutine
58!EOC
59
subroutine ggamt_1(tsh, is, np, rho, grho, g2rho, g3rho)
Definition ggamt_1.f90:10
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