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(:),grfmt(:,:,:),grfir(:,:)
14 complex(8), allocatable :: zvclmt(:,:),zvclir(:)
15 allocate(rfmt(npmtmax,natmtot),rfir(ngtot))
16 allocate(grfmt(npmtmax,natmtot,3),grfir(ngtot,3))
17 allocate(zvclmt(npmtmax,natmtot),zvclir(ngtot))
18 ! compute the divergence of B_xc
19 rfmt(1:npmtmax,1:natmtot)=0.d0
20 rfir(1:ngtot)=0.d0
21 do idm=1,3
22  call gradrf(bxcmt(:,:,idm),bxcir(:,idm),grfmt,grfir)
23  do ias=1,natmtot
24  is=idxis(ias)
25  np=npmt(is)
26  rfmt(1:np,ias)=rfmt(1:np,ias)+grfmt(1:np,ias,idm)
27  end do
28  rfir(1:ngtot)=rfir(1:ngtot)+grfir(1:ngtot,idm)
29 end do
30 ! convert real muffin-tin divergence to complex spherical harmonic expansion
31 do ias=1,natmtot
32  is=idxis(ias)
33  call rtozfmt(nrmt(is),nrmti(is),rfmt(:,ias),zvclmt(:,ias))
34 end do
35 ! store real interstitial divergence in a complex array
36 zvclir(1:ngtot)=rfir(1:ngtot)
37 ! solve the complex Poisson's equation
40  jlgrmt,ylmg,sfacg,npmtmax,zvclmt,zvclir)
41 ! convert complex muffin-tin potential to real spherical harmonic expansion
42 do ias=1,natmtot
43  is=idxis(ias)
44  call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),rfmt(:,ias))
45 end do
46 ! store complex interstitial potential in real array
47 rfir(1:ngtot)=dble(zvclir(1:ngtot))
48 ! compute the gradient
49 call gradrf(rfmt,rfir,grfmt,grfir)
50 ! add gradient over 4π to existing B_xc
51 t1=1.d0/fourpi
52 do idm=1,3
53  do ias=1,natmtot
54  is=idxis(ias)
55  np=npmt(is)
56  bxcmt(1:np,ias,idm)=bxcmt(1:np,ias,idm)+t1*grfmt(1:np,ias,idm)
57  end do
58 end do
59 bxcir(1:ngtot,1:3)=bxcir(1:ngtot,1:3)+t1*grfir(1:ngtot,1:3)
60 deallocate(rfmt,rfir,grfmt,grfir)
61 deallocate(zvclmt,zvclir)
62 end subroutine
63 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
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
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zvclmt)
Definition: genzvclmt.f90:7
integer ngvec
Definition: modmain.f90:396
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
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:1228
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