The Elk Code
xc_c_tb09.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 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 subroutine xc_c_tb09
7 use modmain
8 use libxcifc
9 implicit none
10 ! local variables
11 integer is,ias,i
12 integer nr,nri,ir
13 real(8), parameter :: alpha=-0.012d0, beta=1.023d0
14 real(8) t1
15 ! allocatable arrays
16 real(8), allocatable :: grfmt(:,:,:),grfir(:,:)
17 real(8), allocatable :: rfmt(:,:),rfir(:)
18 ! external functions
19 real(8), external :: rfint
20 ! compute the gradient of the density
21 allocate(grfmt(npmtmax,natmtot,3),grfir(ngtot,3))
22 call gradrf(rhomt,rhoir,grfmt,grfir)
23 allocate(rfmt(npmtmax,natmtot))
24 do ias=1,natmtot
25  is=idxis(ias)
26  nr=nrmt(is)
27  nri=nrmti(is)
28 ! convert muffin-tin density to spherical coordinates
29  call rbsht(nr,nri,rhomt(:,ias),rfmt(:,ias))
30 ! convert muffin-tin gradient to spherical coordinates
31  do i=1,3
32  call rbshtip(nr,nri,grfmt(:,ias,i))
33  end do
34 ! integrand in muffin-tin
35  do i=1,npmt(is)
36  t1=sqrt(grfmt(i,ias,1)**2+grfmt(i,ias,2)**2+grfmt(i,ias,3)**2)
37  rfmt(i,ias)=t1/rfmt(i,ias)
38  end do
39 ! convert to spherical harmonics
40  call rfshtip(nr,nri,rfmt(:,ias))
41 end do
42 deallocate(grfmt)
43 ! integrand in interstitial
44 allocate(rfir(ngtot))
45 do ir=1,ngtot
46  t1=sqrt(grfir(ir,1)**2+grfir(ir,2)**2+grfir(ir,3)**2)
47  rfir(ir)=t1/rhoir(ir)
48 end do
49 ! integrate over the unit cell
50 t1=rfint(rfmt,rfir)
51 ! set the constant
52 c_tb09=alpha+beta*sqrt(abs(t1)/omega)
53 deallocate(grfir,rfmt,rfir)
54 end subroutine
55 
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
integer ngtot
Definition: modmain.f90:390
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine xc_c_tb09
Definition: xc_c_tb09.f90:7
real(8) omega
Definition: modmain.f90:20
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition: gradrf.f90:7
real(8) c_tb09
Definition: modmain.f90:678
subroutine rfshtip(nr, nri, rfmt)
Definition: rfshtip.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer npmtmax
Definition: modmain.f90:216
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150