The Elk Code
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 
6 subroutine dpotcoul
7 use modmain
8 use modphonon
9 implicit none
10 ! local variables
11 integer np
12 ! allocatable arrays
13 complex(8), allocatable :: gzfmt(:,:)
14 ! solve the complex Poisson's equation in the muffin-tins
16 ! calculate the gradient of the nuclear potential
17 allocate(gzfmt(npmtmax,3))
18 call gradzvcln(isph,gzfmt)
19 ! subtract gradient component corresponding to the phonon polarisation
20 np=npmt(isph)
21 dvclmt(1:np,iasph)=dvclmt(1:np,iasph)-gzfmt(1:np,ipph)
22 deallocate(gzfmt)
23 ! solve Poisson's equation in the entire unit cell
26 end subroutine
27 
integer, dimension(3) ngridg
Definition: modmain.f90:386
real(8), dimension(:,:,:), allocatable jlgqrmt
Definition: modphonon.f90:68
complex(8), dimension(:,:), pointer, contiguous drhomt
Definition: modphonon.f90:100
complex(8), dimension(:,:), allocatable dvclmt
Definition: modphonon.f90:104
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
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine dpotcoul
Definition: dpotcoul.f90:7
integer iasph
Definition: modphonon.f90:15
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine gradzvcln(is, gzfmt)
Definition: gradzvcln.f90:7
integer ipph
Definition: modphonon.f90:15
integer ngvec
Definition: modmain.f90:396
integer isph
Definition: modphonon.f90:15
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition: genzvclmt.f90:7
complex(8), dimension(:), pointer, contiguous drhoir
Definition: modphonon.f90:100
complex(8), dimension(:,:), allocatable sfacgq
Definition: modphonon.f90:72
integer npmtmax
Definition: modmain.f90:216
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
integer nrmtmax
Definition: modmain.f90:152
real(8), dimension(:), allocatable gqc
Definition: modphonon.f90:64
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
complex(8), dimension(:), allocatable dvclir
Definition: modphonon.f90:104
real(8), dimension(:), allocatable gclgq
Definition: modphonon.f90:66
complex(8), dimension(:,:), allocatable ylmgq
Definition: modphonon.f90:70
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150