The Elk Code
 
Loading...
Searching...
No Matches
ggamt_2a.f90
Go to the documentation of this file.
1
2! Copyright (C) 2009 T. McQueen and J. K. Dewhurst.
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_2a
8! !INTERFACE:
9subroutine ggamt_2a(tsh,is,np,rho,g2rho,gvrho,grho2)
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Spin-unpolarised version of {\tt ggamt\_sp\_2a}.
14!
15! !REVISION HISTORY:
16! Created November 2009 (JKD and TMcQ)
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) :: g2rho(np),gvrho(np,3),grho2(np)
25! local variables
26integer nr,nri,i
27! automatic arrays
28real(8) grfmt(np,3),rfmt1(np),rfmt2(np)
29nr=nrmt(is)
30nri=nrmti(is)
31! compute grad^2 rho in spherical coordinates
32if (tsh) then
33 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rho,rfmt2)
34else
35 call rfsht(nr,nri,rho,rfmt1)
36 call grad2rfmt(nr,nri,rlmt(:,-1,is),rlmt(:,-2,is),wcrmt(:,:,is),rfmt1,rfmt2)
37end if
38call rbsht(nr,nri,rfmt2,g2rho)
39! compute grad rho in spherical coordinates
40if (tsh) then
41 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rho,np,grfmt)
42else
43 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt1,np,grfmt)
44end if
45do i=1,3
46 call rbsht(nr,nri,grfmt(:,i),gvrho(:,i))
47end do
48! (grad rho)^2
49grho2(1:np)=gvrho(1:np,1)**2+gvrho(1:np,2)**2+gvrho(1:np,3)**2
50end subroutine
51!EOC
52
subroutine ggamt_2a(tsh, is, np, rho, g2rho, gvrho, grho2)
Definition ggamt_2a.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