The Elk Code
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 
6 subroutine gengvnsmt
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer is,ias,l
12 integer nr,nri,iro
13 integer np,i0,i1
14 ! automatic arrays
15 complex(8) zvclmt(npmtmax),gzfmt(npmtmax,3)
16 ! allocate global array
17 if (allocated(gvnsmt)) deallocate(gvnsmt)
18 allocate(gvnsmt(npmtmax,3,natmtot))
19 ! loop over atoms
20 do ias=1,natmtot
21  is=idxis(ias)
22  nr=nrmt(is)
23  nri=nrmti(is)
24  iro=nri+1
25  np=npmt(is)
26  do l=1,3
27 ! convert static density to complex spherical harmonic expansion
28  call rtozfmt(nr,nri,rhosmt(:,ias,l),zvclmt)
29 ! solve Poisson's equation in the muffin-tin
30  call zpotclmt(nr,nri,nrmtmax,rlmt(:,:,is),wprmt(:,:,is),zvclmt)
31 ! add the nuclear Coulomb potential
32  i1=lmmaxi*(nri-1)+1
33  zvclmt(1:i1:lmmaxi)=zvclmt(1:i1:lmmaxi)+vcln(1:nri,is)
34  i0=i1+lmmaxi
35  i1=lmmaxo*(nr-iro)+i0
36  zvclmt(i0:i1:lmmaxo)=zvclmt(i0:i1:lmmaxo)+vcln(iro:nr,is)
37 ! compute the gradient of the potential and store in global array
38  call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),zvclmt,npmtmax,gzfmt)
39  gvnsmt(1:np,l,ias)=gzfmt(1:np,l)
40  end do
41 end do
42 end subroutine
43 
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
complex(8), dimension(:,:,:), allocatable gvnsmt
Definition: modtddft.f90:88
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
pure subroutine zpotclmt(nr, nri, ld, rl, wpr, zvclmt)
Definition: zpotclmt.f90:10
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
subroutine gengvnsmt
Definition: gengvnsmt.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), dimension(:,:,:), allocatable rhosmt
Definition: modtddft.f90:82
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150