The Elk Code
 
Loading...
Searching...
No Matches
gengvnsmt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine gengvnsmt
7use modmain
8use modtddft
9implicit none
10! local variables
11integer is,ias
12integer nr,nri,iro
13integer np,i0,i1
14! allocatable arrays
15real(8), allocatable :: rfmt(:)
16complex(8), allocatable :: zrhomt(:),zvclmt(:)
17! allocate local arrays
18allocate(rfmt(npmtmax),zrhomt(npmtmax),zvclmt(npmtmax))
19! allocate global array
20if (allocated(gvnsmt)) deallocate(gvnsmt)
21allocate(gvnsmt(npmtmax,3,natmtot))
22! loop over atoms
23do ias=1,natmtot
24 is=idxis(ias)
25 nr=nrmt(is)
26 nri=nrmti(is)
27 iro=nri+1
28 np=npmt(is)
29! take the average of the static density over the three directions
30 rfmt(1:np)=(1.d0/3.d0) &
31 *(rhosmt(1:np,ias,1)+rhosmt(1:np,ias,2)+rhosmt(1:np,ias,3))
32! convert to complex spherical harmonics expansion
33 call rtozfmt(nr,nri,rfmt,zrhomt)
34! solve Poisson's equation in the muffin-tin
35 call zpotclmt(nr,nri,nrmtmax,rlmt(:,:,is),wprmt(:,:,is),zrhomt,zvclmt)
36! add the nuclear Coulomb potential
37 i1=lmmaxi*(nri-1)+1
38 zvclmt(1:i1:lmmaxi)=zvclmt(1:i1:lmmaxi)+vcln(1:nri,is)
39 i0=i1+lmmaxi
40 i1=lmmaxo*(nr-iro)+i0
41 zvclmt(i0:i1:lmmaxo)=zvclmt(i0:i1:lmmaxo)+vcln(iro:nr,is)
42! compute the gradient of the potential and store in global array
43 call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),zvclmt,npmtmax, &
44 gvnsmt(:,:,ias))
45end do
46deallocate(rfmt,zrhomt,zvclmt)
47end subroutine
48
subroutine gengvnsmt
Definition gengvnsmt.f90:7
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.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
integer lmmaxi
Definition modmain.f90:207
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:), allocatable wprmt
Definition modmain.f90:185
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vcln
Definition modmain.f90:97
integer nrmtmax
Definition modmain.f90:152
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
real(8), dimension(:,:,:), allocatable rhosmt
Definition modtddft.f90:82
complex(8), dimension(:,:,:), allocatable gvnsmt
Definition modtddft.f90:88
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
pure subroutine zpotclmt(nr, nri, ld, rl, wpr, zrhomt, zvclmt)
Definition zpotclmt.f90:10