The Elk Code
projsbf.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 projsbf
7 use modmain
8 implicit none
9 ! local variables
10 integer idm,is,ias,np
11 real(8) t1
12 ! allocatable arrays
13 real(8), allocatable :: rfmt(:,:),rfir(:)
14 real(8), allocatable :: grfmt(:,:,:),grfir(:,:)
15 complex(8), allocatable :: zrhomt(:,:),zrhoir(:)
16 complex(8), allocatable :: zvclmt(:,:),zvclir(:)
17 allocate(rfmt(npmtmax,natmtot),rfir(ngtot))
18 allocate(grfmt(npmtmax,natmtot,3),grfir(ngtot,3))
19 allocate(zrhomt(npmtmax,natmtot),zrhoir(ngtot))
20 allocate(zvclmt(npmtmax,natmtot),zvclir(ngtot))
21 ! compute the divergence of B_xc
22 rfmt(1:npmtmax,1:natmtot)=0.d0
23 rfir(1:ngtot)=0.d0
24 do idm=1,3
25  call gradrf(bxcmt(:,:,idm),bxcir(:,idm),grfmt,grfir)
26  do ias=1,natmtot
27  is=idxis(ias)
28  np=npmt(is)
29  rfmt(1:np,ias)=rfmt(1:np,ias)+grfmt(1:np,ias,idm)
30  end do
31  rfir(1:ngtot)=rfir(1:ngtot)+grfir(1:ngtot,idm)
32 end do
33 ! convert real muffin-tin divergence to complex spherical harmonic expansion
34 do ias=1,natmtot
35  is=idxis(ias)
36  call rtozfmt(nrmt(is),nrmti(is),rfmt(:,ias),zrhomt(:,ias))
37 end do
38 ! store real interstitial divergence in a complex array
39 zrhoir(1:ngtot)=rfir(1:ngtot)
40 ! solve the complex Poisson's equation
41 call genzvclmt(nrmt,nrmti,nrmtmax,rlmt,wprmt,npmtmax,zrhomt,zvclmt)
43  jlgrmt,ylmg,sfacg,zrhoir,npmtmax,zvclmt,zvclir)
44 ! convert complex muffin-tin potential to real spherical harmonic expansion
45 do ias=1,natmtot
46  is=idxis(ias)
47  call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),rfmt(:,ias))
48 end do
49 ! store complex interstitial potential in real array
50 rfir(1:ngtot)=dble(zvclir(1:ngtot))
51 ! compute the gradient
52 call gradrf(rfmt,rfir,grfmt,grfir)
53 ! add gradient over 4π to existing B_xc
54 t1=1.d0/fourpi
55 do idm=1,3
56  do ias=1,natmtot
57  is=idxis(ias)
58  np=npmt(is)
59  bxcmt(1:np,ias,idm)=bxcmt(1:np,ias,idm)+t1*grfmt(1:np,ias,idm)
60  end do
61 end do
62 bxcir(1:ngtot,1:3)=bxcir(1:ngtot,1:3)+t1*grfir(1:ngtot,1:3)
63 deallocate(rfmt,rfir,grfmt,grfir)
64 deallocate(zrhomt,zrhoir,zvclmt,zvclir)
65 end subroutine
66 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, zrhoir, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition: gradrf.f90:7
complex(8), dimension(:,:), allocatable ylmg
Definition: modmain.f90:428
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
real(8), dimension(:,:,:), allocatable bxcmt
Definition: modmain.f90:636
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition: genzvclmt.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition: ztorfmt.f90:7
real(8), dimension(:), allocatable gclg
Definition: modmain.f90:424
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
integer npmtmax
Definition: modmain.f90:216
real(8), parameter fourpi
Definition: modmain.f90:1234
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
subroutine projsbf
Definition: projsbf.f90:7
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150