The Elk Code
 
Loading...
Searching...
No Matches
potcoul.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: potcoul
8! !INTERFACE:
9subroutine potcoul
10! !USES:
11use modmain
12use modomp
13! !DESCRIPTION:
14! Calculates the Coulomb potential of the real charge density stored in the
15! global variables {\tt rhomt} and {\tt rhoir} by solving Poisson's equation.
16! These variables are coverted to complex representations and passed to the
17! routine {\tt zpotcoul}.
18!
19! !REVISION HISTORY:
20! Created April 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! local variables
25integer is,ias,nthd
26integer nr,nri,iro,i0,i1
27! allocatable arrays
28complex(8), allocatable :: zrhomt(:,:),zrhoir(:)
29complex(8), allocatable :: zvclmt(:,:),zvclir(:)
30allocate(zrhomt(npmtmax,natmtot))
31! convert real muffin-tin charge density to complex spherical harmonic expansion
32call holdthd(natmtot,nthd)
33!$OMP PARALLEL DO DEFAULT(SHARED) &
34!$OMP PRIVATE(is) &
35!$OMP SCHEDULE(DYNAMIC) &
36!$OMP NUM_THREADS(nthd)
37do ias=1,natmtot
38 is=idxis(ias)
39 call rtozfmt(nrmt(is),nrmti(is),rhomt(:,ias),zrhomt(:,ias))
40end do
41!$OMP END PARALLEL DO
42call freethd(nthd)
43! solve the complex Poisson's equation in the muffin-tins
44allocate(zvclmt(npmtmax,natmtot))
45call genzvclmt(nrmt,nrmti,nrmtmax,rlmt,wprmt,npmtmax,zrhomt,zvclmt)
46deallocate(zrhomt)
47! add the nuclear monopole potentials
48do ias=1,natmtot
49 is=idxis(ias)
50 nr=nrmt(is)
51 nri=nrmti(is)
52 iro=nri+1
53 i1=lmmaxi*(nri-1)+1
54 zvclmt(1:i1:lmmaxi,ias)=zvclmt(1:i1:lmmaxi,ias)+vcln(1:nri,is)
55 i0=i1+lmmaxi
56 i1=lmmaxo*(nr-iro)+i0
57 zvclmt(i0:i1:lmmaxo,ias)=zvclmt(i0:i1:lmmaxo,ias)+vcln(iro:nr,is)
58end do
59! apply atomic displacement potential if required
60if (tatdisp) call zvcldisp(zvclmt)
61allocate(zrhoir(ngtot),zvclir(ngtot))
62! store real interstitial charge density in complex array
63zrhoir(1:ngtot)=rhoir(1:ngtot)
64! solve Poisson's equation in the entire unit cell
66 jlgrmt,ylmg,sfacg,zrhoir,npmtmax,zvclmt,zvclir)
67deallocate(zrhoir)
68! convert complex muffin-tin potential to real spherical harmonic expansion
69call holdthd(natmtot+1,nthd)
70!$OMP PARALLEL DEFAULT(SHARED) &
71!$OMP PRIVATE(ias,is) &
72!$OMP NUM_THREADS(nthd)
73!$OMP DO SCHEDULE(DYNAMIC)
74do ias=1,natmtot
75 is=idxis(ias)
76 call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),vclmt(:,ias))
77end do
78!$OMP END DO NOWAIT
79! store complex interstitial potential in real array
80!$OMP SINGLE
81vclir(1:ngtot)=dble(zvclir(1:ngtot))
82!$OMP END SINGLE
83!$OMP END PARALLEL
84call freethd(nthd)
85deallocate(zvclmt,zvclir)
86! apply constant electric field if required
87if (tefield) call potefield
88end subroutine
89!EOC
90
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition genzvclmt.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
real(8), dimension(:,:,:), allocatable jlgrmt
Definition modmain.f90:426
integer lmmaxi
Definition modmain.f90:207
integer ngvec
Definition modmain.f90:396
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
complex(8), dimension(:,:), allocatable sfacg
Definition modmain.f90:430
logical tefield
Definition modmain.f90:310
integer natmtot
Definition modmain.f90:40
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:), allocatable vclir
Definition modmain.f90:624
logical tatdisp
Definition modmain.f90:59
real(8), dimension(:,:,:), allocatable wprmt
Definition modmain.f90:185
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition modmain.f90:624
real(8), dimension(:,:), allocatable vcln
Definition modmain.f90:97
real(8), dimension(:), allocatable gclg
Definition modmain.f90:424
integer nrmtmax
Definition modmain.f90:152
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
complex(8), dimension(:,:), allocatable ylmg
Definition modmain.f90:428
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine potcoul
Definition potcoul.f90:10
subroutine potefield
Definition potefield.f90:7
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
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
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition ztorfmt.f90:7
subroutine zvcldisp(zvclmt)
Definition zvcldisp.f90:7