forked from qsnake/elk
/
potcoul.f90
87 lines (84 loc) · 2.49 KB
/
potcoul.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.
!BOP
! !ROUTINE: potcoul
! !INTERFACE:
subroutine potcoul
! !USES:
use modmain
! !DESCRIPTION:
! Calculates the Coulomb potential of the real charge density stored in the
! global variables {\tt rhomt} and {\tt rhoir} by solving Poisson's equation.
! These variables are coverted to complex representations and passed to the
! routine {\tt zpotcoul}.
!
! !REVISION HISTORY:
! Created April 2003 (JKD)
!EOP
!BOC
implicit none
! local variables
integer is,ia,ias,ir,lmax
real(8) t1
complex(8) zrho0
! automatic arrays
real(8) vn(nrmtmax)
! allocatable arrays
real(8), allocatable :: jlgr(:,:,:)
complex(8), allocatable :: zrhomt(:,:,:)
complex(8), allocatable :: zrhoir(:)
complex(8), allocatable :: zvclmt(:,:,:)
complex(8), allocatable :: zvclir(:)
allocate(jlgr(0:lmaxvr+npsden+1,ngvec,nspecies))
allocate(zrhomt(lmmaxvr,nrmtmax,natmtot))
allocate(zrhoir(ngrtot))
allocate(zvclmt(lmmaxvr,nrmtmax,natmtot))
allocate(zvclir(ngrtot))
! convert real muffin-tin charge density to complex spherical harmonic expansion
do is=1,nspecies
do ia=1,natoms(is)
ias=idxas(ia,is)
do ir=1,nrmt(is)
call rtozflm(lmaxvr,rhomt(:,ir,ias),zrhomt(:,ir,ias))
end do
end do
end do
! store real interstitial charge density in complex array
zrhoir(:)=rhoir(:)
! compute the required spherical Bessel functions
lmax=lmaxvr+npsden+1
call genjlgpr(lmax,gc,jlgr)
! solve the complex Poisson's equation in the muffin-tins
call genzvclmt(nrmt,spnrmax,spr,nrmtmax,zrhomt,zvclmt)
! add the nuclear monopole potentials
t1=1.d0/y00
do is=1,nspecies
call potnucl(ptnucl,nrmt(is),spr(:,is),spzn(is),vn)
do ia=1,natoms(is)
ias=idxas(ia,is)
do ir=1,nrmt(is)
zvclmt(1,ir,ias)=zvclmt(1,ir,ias)+t1*vn(ir)
end do
end do
end do
! solve Poisson's equation in the entire unit cell
call zpotcoul(nrmt,spnrmax,spr,1,gc,jlgr,ylmg,sfacg,zrhoir,nrmtmax,zvclmt, &
zvclir,zrho0)
! convert complex muffin-tin potential to real spherical harmonic expansion
do is=1,nspecies
do ia=1,natoms(is)
ias=idxas(ia,is)
do ir=1,nrmt(is)
call ztorflm(lmaxvr,zvclmt(:,ir,ias),vclmt(:,ir,ias))
end do
end do
end do
! store complex interstitial potential in real array
vclir(:)=dble(zvclir(:))
! apply constant electric field if required
if (efieldpol) call potefield
deallocate(jlgr,zrhomt,zrhoir,zvclmt,zvclir)
return
end subroutine
!EOC