The Elk Code
potcoulu.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 T. Mueller, 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 potcoulu
7 use modmain
8 use modulr
9 use modmpi
10 use modomp
11 implicit none
12 ! local variables
13 integer ifq,is,ias
14 integer nr,nri,ir
15 integer nrc,nrci,i
16 integer n,lp,nthd
17 ! allocatable arrays
18 complex(8), allocatable :: zvclmt(:,:)
19 call holdthd(nfqrz/np_mpi,nthd)
20 !$OMP PARALLEL DEFAULT(SHARED) &
21 !$OMP PRIVATE(zvclmt,ias,is,nr,nri) &
22 !$OMP PRIVATE(i,ir,nrc,nrci) &
23 !$OMP NUM_THREADS(nthd)
24 allocate(zvclmt(npmtmax,natmtot))
25 !$OMP DO SCHEDULE(DYNAMIC)
26 do ifq=1,nfqrz
27 ! distribute among MPI processes
28  if (mod(ifq-1,np_mpi) /= lp_mpi) cycle
29 ! convert the complex muffin-tin density from coarse to fine radial mesh
30  do ias=1,natmtot
31  is=idxis(ias)
32  zvclmt(1:npcmt(is),ias)=rhoqmt(1:npcmt(is),ias,ifq)
33  end do
34  call zfmtctof(zvclmt)
35 ! solve the complex Poisson's equation in the muffin-tins
37 ! add the nuclear monopole potentials for Q=0
38  if (ifq == 1) then
39  do ias=1,natmtot
40  is=idxis(ias)
41  nr=nrmt(is)
42  nri=nrmti(is)
43  i=1
44  do ir=1,nri
45  zvclmt(i,ias)=zvclmt(i,ias)+vcln(ir,is)
46  i=i+lmmaxi
47  end do
48  do ir=nri+1,nr
49  zvclmt(i,ias)=zvclmt(i,ias)+vcln(ir,is)
50  i=i+lmmaxo
51  end do
52  end do
53  end if
54 ! convert the interstitial density from coarse to fine grid
55  call zfirctof(rhoqir(:,ifq),vsqir(:,ifq))
56 ! solve Poisson's equation in the entire unit cell
58  gclgq(:,ifq),ngvec,jlgqrmt(:,:,:,ifq),ylmgq(:,:,ifq),sfacgq(:,:,ifq), &
59  npmtmax,zvclmt,vsqir(:,ifq))
60  do ias=1,natmtot
61  is=idxis(ias)
62  nrc=nrcmt(is)
63  nrci=nrcmti(is)
64 ! convert from fine to coarse radial mesh
65  call zfmtftoc(nrc,nrci,zvclmt(:,ias),vsqmt(:,ias,ifq))
66 ! convert to spherical coordinates
67  call zbshtip(nrc,nrci,vsqmt(:,ias,ifq))
68  end do
69 end do
70 !$OMP END DO
71 deallocate(zvclmt)
72 !$OMP END PARALLEL
73 call freethd(nthd)
74 ! broadcast potentials to every MPI process
75 if (np_mpi > 1) then
77  do ifq=1,nfqrz
78  lp=mod(ifq-1,np_mpi)
79  call mpi_bcast(vsqmt(:,:,ifq),n,mpi_double_complex,lp,mpicom,ierror)
80  end do
81  do ifq=1,nfqrz
82  lp=mod(ifq-1,np_mpi)
83  call mpi_bcast(vsqir(:,ifq),ngtot,mpi_double_complex,lp,mpicom,ierror)
84  end do
85 end if
86 end subroutine
87 
subroutine zfmtctof(zfmt)
Definition: zfmtctof.f90:7
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition: modulr.f90:82
complex(8), dimension(:,:,:), allocatable sfacgq
Definition: modulr.f90:44
real(8), dimension(:,:,:,:), allocatable jlgqrmt
Definition: modulr.f90:40
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer npcmtmax
Definition: modmain.f90:216
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
integer lmmaxo
Definition: modmain.f90:203
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition: modulr.f90:82
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
Definition: modomp.f90:6
integer np_mpi
Definition: modmpi.f90:13
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine potcoulu
Definition: potcoulu.f90:7
subroutine zbshtip(nr, nri, zfmt)
Definition: zbshtip.f90:7
pure subroutine zfmtftoc(nrc, nrci, zfmt, zfcmt)
Definition: zfmtftoc.f90:7
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
complex(8), dimension(:,:,:), allocatable ylmgq
Definition: modulr.f90:42
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:), allocatable gqc
Definition: modulr.f90:36
integer lmmaxi
Definition: modmain.f90:207
integer nfqrz
Definition: modmain.f90:539
Definition: modmpi.f90:6
subroutine zfirctof(zfirc, zfir)
Definition: zfirctof.f90:7
integer lp_mpi
Definition: modmpi.f90:15
complex(8), dimension(:,:,:), allocatable rhoqmt
Definition: modulr.f90:58
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
real(8), dimension(:,:), allocatable gclgq
Definition: modulr.f90:38
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
Definition: modulr.f90:6
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
complex(8), dimension(:,:), allocatable rhoqir
Definition: modulr.f90:58