The Elk Code
 
Loading...
Searching...
No Matches
gradzvcln.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 gradzvcln(is,gzfmt)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: is
11complex(8), intent(out) :: gzfmt(npmtmax,3)
12! local variables
13integer nr,nri,iro,i0,i1
14! allocatable arrays
15complex(8), allocatable :: zvclmt(:)
16allocate(zvclmt(npmtmax))
17nr=nrmt(is)
18nri=nrmti(is)
19iro=nri+1
20! convert nuclear Coulomb potential to complex spherical harmonics expansion
21zvclmt(1:npmt(is))=0.d0
22i1=lmmaxi*(nri-1)+1
23zvclmt(1:i1:lmmaxi)=vcln(1:nri,is)
24i0=i1+lmmaxi
25i1=lmmaxo*(nr-iro)+i0
26zvclmt(i0:i1:lmmaxo)=vcln(iro:nr,is)
27! compute the gradient of the potential
28call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),zvclmt,npmtmax,gzfmt)
29deallocate(zvclmt)
30end subroutine
31
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
subroutine gradzvcln(is, gzfmt)
Definition gradzvcln.f90:7
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(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vcln
Definition modmain.f90:97
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179