The Elk Code
 
Loading...
Searching...
No Matches
dpotcoul.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine dpotcoul
7use modmain
8use modphonon
9implicit none
10! local variables
11integer np
12! allocatable arrays
13complex(8), allocatable :: gzfmt(:,:)
14! solve the complex Poisson's equation in the muffin-tins
16! calculate the gradient of the nuclear potential
17allocate(gzfmt(npmtmax,3))
18call gradzvcln(isph,gzfmt)
19! subtract gradient component corresponding to the phonon polarisation
20np=npmt(isph)
21dvclmt(1:np,iasph)=dvclmt(1:np,iasph)-gzfmt(1:np,ipph)
22deallocate(gzfmt)
23! solve Poisson's equation in the entire unit cell
26end subroutine
27
subroutine dpotcoul
Definition dpotcoul.f90:7
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition genzvclmt.f90:7
subroutine gradzvcln(is, gzfmt)
Definition gradzvcln.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer, dimension(3) ngridg
Definition modmain.f90:386
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer ngvec
Definition modmain.f90:396
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:,:,:), allocatable wprmt
Definition modmain.f90:185
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
integer nrmtmax
Definition modmain.f90:152
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
complex(8), dimension(:,:), allocatable ylmgq
Definition modphonon.f90:70
real(8), dimension(:), allocatable gclgq
Definition modphonon.f90:66
integer ipph
Definition modphonon.f90:15
real(8), dimension(:,:,:), allocatable jlgqrmt
Definition modphonon.f90:68
real(8), dimension(:), allocatable gqc
Definition modphonon.f90:64
integer iasph
Definition modphonon.f90:15
integer isph
Definition modphonon.f90:15
complex(8), dimension(:), allocatable drhoir
Definition modphonon.f90:98
complex(8), dimension(:,:), allocatable drhomt
Definition modphonon.f90:98
complex(8), dimension(:,:), allocatable sfacgq
Definition modphonon.f90:72
complex(8), dimension(:,:), allocatable dvclmt
complex(8), dimension(:), allocatable dvclir
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