The Elk Code
 
Loading...
Searching...
No Matches
genfxcr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 S. Sharma, J. K. Dewhurst 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
6subroutine genfxcr(tsh,fxcmt,fxcir)
7use modmain
8use modtddft
9use modfxcifc
10implicit none
11! arguments
12logical, intent(in) :: tsh
13real(8), intent(out) :: fxcmt(npmtmax,natmtot),fxcir(ngtot)
14! local variables
15integer idm,is,ias
16integer nr,nri,ir,np,i,n
17real(8) t1
18real(8), allocatable :: rho(:),rhoup(:),rhodn(:),mag(:,:)
19real(8), allocatable :: fxc(:),fxcuu(:),fxcud(:),fxcdd(:)
20! number of independent spin components
21n=npmtmax
22allocate(rho(n),fxc(n))
23if (spinpol) then
24 allocate(mag(n,3))
25 n=max(n,ngtot)
26 allocate(rhoup(n),rhodn(n))
27 allocate(fxcuu(n),fxcud(n),fxcdd(n))
28end if
29!---------------------------!
30! muffin-tin kernel !
31!---------------------------!
32do ias=1,natmtot
33 is=idxis(ias)
34 nr=nrmt(is)
35 nri=nrmti(is)
36 np=npmt(is)
37! compute the density in spherical coordinates
38 call rbsht(nr,nri,rhomt(:,ias),rho)
39 if (spinpol) then
40!------------------------!
41! spin-polarised !
42!------------------------!
43! magnetisation in spherical coordinates
44 do idm=1,ndmag
45 call rbsht(nr,nri,magmt(:,ias,idm),mag(:,idm))
46 end do
47 if (ncmag) then
48! non-collinear (use Kubler's trick)
49 do i=1,np
50! compute rhoup=(rho+|m|)/2 and rhodn=(rho-|m|)/2
51 t1=sqrt(mag(i,1)**2+mag(i,2)**2+mag(i,3)**2)
52 rhoup(i)=0.5d0*(rho(i)+t1)
53 rhodn(i)=0.5d0*(rho(i)-t1)
54 end do
55 else
56! collinear
57 do i=1,np
58! compute rhoup=(rho+m_z)/2 and rhodn=(rho-m_z)/2
59 rhoup(i)=0.5d0*(rho(i)+mag(i,1))
60 rhodn(i)=0.5d0*(rho(i)-mag(i,1))
61 end do
62 end if
63! compute fxc
64 call fxcifc(fxctype,n=np,rhoup=rhoup,rhodn=rhodn,fxcuu=fxcuu,fxcud=fxcud, &
65 fxcdd=fxcdd)
66! form the scalar quantity dv/drho
67 do i=1,np
68 fxc(i)=0.25d0*(fxcuu(i)+2.d0*fxcud(i)+fxcdd(i))
69 end do
70 else
71!--------------------------!
72! spin-unpolarised !
73!--------------------------!
74 call fxcifc(fxctype,n=np,rho=rho,fxc=fxc)
75 end if
76 if (tsh) then
77! convert fxc to spherical harmonics if required
78 call rfsht(nr,nri,fxc,fxcmt(:,ias))
79 else
80 fxcmt(1:np,ias)=fxc(1:np)
81 end if
82end do
83!-----------------------------!
84! interstitial kernel !
85!-----------------------------!
86if (spinpol) then
87!------------------------!
88! spin-polarised !
89!------------------------!
90 if (ncmag) then
91! non-collinear
92 do ir=1,ngtot
93 t1=sqrt(magir(ir,1)**2+magir(ir,2)**2+magir(ir,3)**2)
94 rhoup(ir)=0.5d0*(rhoir(ir)+t1)
95 rhodn(ir)=0.5d0*(rhoir(ir)-t1)
96 end do
97 else
98! collinear
99 do ir=1,ngtot
100 rhoup(ir)=0.5d0*(rhoir(ir)+magir(ir,1))
101 rhodn(ir)=0.5d0*(rhoir(ir)-magir(ir,1))
102 end do
103 end if
104! compute fxc
105 call fxcifc(fxctype,n=ngtot,rhoup=rhoup,rhodn=rhodn,fxcuu=fxcuu,fxcud=fxcud, &
106 fxcdd=fxcdd)
107 do ir=1,ngtot
108 fxcir(ir)=0.25d0*(fxcuu(ir)+2.d0*fxcud(ir)+fxcdd(ir))
109 end do
110else
111!--------------------------!
112! spin-unpolarised !
113!--------------------------!
114 call fxcifc(fxctype,n=ngtot,rho=rhoir,fxc=fxcir)
115end if
116deallocate(rho,fxc)
117if (spinpol) then
118 deallocate(mag,rhoup,rhodn)
119 deallocate(fxcuu,fxcud,fxcdd)
120end if
121end subroutine
122
subroutine genfxcr(tsh, fxcmt, fxcir)
Definition genfxcr.f90:7
subroutine fxcifc(fxctype, n, rho, rhoup, rhodn, fxc, fxcuu, fxcud, fxcdd)
Definition modfxcifc.f90:13
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
logical ncmag
Definition modmain.f90:240
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
logical spinpol
Definition modmain.f90:228
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
integer ndmag
Definition modmain.f90:238
integer, dimension(3) fxctype
Definition modtddft.f90:12
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition rfsht.f90:7