The Elk Code
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:
9 subroutine potcoul
10 ! !USES:
11 use modmain
12 use 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
23 implicit none
24 ! local variables
25 integer is,ias,nthd
26 integer nr,nri,iro,i0,i1
27 ! allocatable arrays
28 complex(8), allocatable :: zvclmt(:,:),zvclir(:)
29 allocate(zvclmt(npmtmax,natmtot))
30 ! convert real muffin-tin charge density to complex spherical harmonic expansion
31 call holdthd(natmtot,nthd)
32 !$OMP PARALLEL DO DEFAULT(SHARED) &
33 !$OMP PRIVATE(is) &
34 !$OMP SCHEDULE(DYNAMIC) &
35 !$OMP NUM_THREADS(nthd)
36 do ias=1,natmtot
37  is=idxis(ias)
38  call rtozfmt(nrmt(is),nrmti(is),rhomt(:,ias),zvclmt(:,ias))
39 end do
40 !$OMP END PARALLEL DO
41 call freethd(nthd)
42 ! solve the complex Poisson's equation in the muffin-tins
44 ! add the nuclear monopole potentials
45 do ias=1,natmtot
46  is=idxis(ias)
47  nr=nrmt(is)
48  nri=nrmti(is)
49  iro=nri+1
50  i1=lmmaxi*(nri-1)+1
51  zvclmt(1:i1:lmmaxi,ias)=zvclmt(1:i1:lmmaxi,ias)+vcln(1:nri,is)
52  i0=i1+lmmaxi
53  i1=lmmaxo*(nr-iro)+i0
54  zvclmt(i0:i1:lmmaxo,ias)=zvclmt(i0:i1:lmmaxo,ias)+vcln(iro:nr,is)
55 end do
56 ! apply atomic displacement potential if required
57 if (tatdisp) call zvcldisp(zvclmt)
58 allocate(zvclir(ngtot))
59 ! store real interstitial charge density in complex array
60 zvclir(1:ngtot)=rhoir(1:ngtot)
61 ! solve Poisson's equation in the entire unit cell
63  jlgrmt,ylmg,sfacg,npmtmax,zvclmt,zvclir)
64 ! convert complex muffin-tin potential to real spherical harmonic expansion
65 call holdthd(natmtot+1,nthd)
66 !$OMP PARALLEL DEFAULT(SHARED) &
67 !$OMP PRIVATE(ias,is) &
68 !$OMP NUM_THREADS(nthd)
69 !$OMP DO SCHEDULE(DYNAMIC)
70 do ias=1,natmtot
71  is=idxis(ias)
72  call ztorfmt(nrmt(is),nrmti(is),zvclmt(:,ias),vclmt(:,ias))
73 end do
74 !$OMP END DO NOWAIT
75 ! store complex interstitial potential in real array
76 !$OMP SINGLE
77 vclir(1:ngtot)=dble(zvclir(1:ngtot))
78 !$OMP END SINGLE
79 !$OMP END PARALLEL
80 call freethd(nthd)
81 deallocate(zvclmt,zvclir)
82 ! apply constant electric field if required
83 if (tefield) call potefield
84 end subroutine
85 !EOC
86 
complex(8), dimension(:,:), allocatable sfacg
Definition: modmain.f90:430
integer, dimension(3) ngridg
Definition: modmain.f90:386
logical tatdisp
Definition: modmain.f90:59
integer ngtot
Definition: modmain.f90:390
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
real(8), dimension(:,:), allocatable vcln
Definition: modmain.f90:97
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
subroutine zvcldisp(zvclmt)
Definition: zvcldisp.f90:7
Definition: modomp.f90:6
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
subroutine potefield
Definition: potefield.f90:7
complex(8), dimension(:,:), allocatable ylmg
Definition: modmain.f90:428
logical tefield
Definition: modmain.f90:310
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zvclmt)
Definition: genzvclmt.f90:7
integer ngvec
Definition: modmain.f90:396
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, ld3, zvclmt, zvclir)
Definition: zpotcoul.f90:11
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:,:), allocatable jlgrmt
Definition: modmain.f90:426
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition: ztorfmt.f90:7
real(8), dimension(:), allocatable gclg
Definition: modmain.f90:424
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer npmtmax
Definition: modmain.f90:216
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine potcoul
Definition: potcoul.f90:10
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150