The Elk Code
 
Loading...
Searching...
No Matches
rhoinit.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: rhoinit
8! !INTERFACE:
9subroutine rhoinit
10! !USES:
11use modmain
12use modomp
13! !DESCRIPTION:
14! Initialises the crystal charge density. Inside the muffin-tins it is set to
15! the spherical atomic density. In the interstitial region it is taken to be
16! constant such that the total charge is correct. Requires that the atomic
17! densities have already been calculated.
18!
19! !REVISION HISTORY:
20! Created January 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! local variables
25integer lmax,is,ia,ias,nthd
26integer nr,nri,nro,nrs,iro,ir
27integer nrc,nrci,irco,irc
28integer l,lm,i0,i1,ig,ifg
29real(8) x,sm,t1,t2
30complex(8) z1,z2
31! automatic arrays
32real(8) ffg(ngvc),wr(nrspmax),jl(0:lmaxi,nrcmtmax)
33complex(8) zfmt(npcmtmax),zfft(ngtc)
34lmax=min(lmaxi,1)
35! compute the superposition of all the atomic density tails
36zfft(1:ngtc)=0.d0
37call holdthd(nspecies,nthd)
38!$OMP PARALLEL DO DEFAULT(SHARED) &
39!$OMP PRIVATE(ffg,wr,nr,nrs,nro,ig) &
40!$OMP PRIVATE(t1,t2,sm,ir,x,ia,ias,ifg) &
41!$OMP REDUCTION(+:zfft) &
42!$OMP SCHEDULE(DYNAMIC) &
43!$OMP NUM_THREADS(nthd)
44do is=1,nspecies
45 nr=nrmt(is)
46 nrs=nrsp(is)
47 nro=nrs-nr+1
48! determine the weights for the radial integral
49 call wsplint(nro,rsp(nr,is),wr(nr))
50 do ig=1,ngvc
51 t1=gc(ig)
52! spherical bessel function j_0(x) times the atomic density tail
53 if (t1 > epslat) then
54 t2=1.d0/t1
55 sm=0.d0
56 do ir=nr,nrs
57 x=t1*rsp(ir,is)
58 sm=sm+t2*sin(x)*rhosp(ir,is)*rsp(ir,is)*wr(ir)
59 end do
60 else
61 sm=sum(rhosp(nr:nrs,is)*(rsp(nr:nrs,is)**2)*wr(nr:nrs))
62 end if
63! apply low-pass filter
64 t1=sm*exp(-4.d0*(gc(ig)/gmaxvr)**2)
65 ffg(ig)=(fourpi/omega)*t1
66 end do
67 do ia=1,natoms(is)
68 ias=idxas(ia,is)
69 do ig=1,ngvc
70 ifg=igfc(ig)
71 zfft(ifg)=zfft(ifg)+ffg(ig)*conjg(sfacg(ig,ias))
72 end do
73 end do
74end do
75!$OMP END PARALLEL DO
76call freethd(nthd)
77! compute the tails in each muffin-tin
78call holdthd(natmtot,nthd)
79!$OMP PARALLEL DO DEFAULT(SHARED) &
80!$OMP PRIVATE(jl,zfmt,is,nrc,nrci) &
81!$OMP PRIVATE(irco,ig,ifg,irc,x) &
82!$OMP PRIVATE(z1,z2,lm,l,i0,i1) &
83!$OMP SCHEDULE(DYNAMIC) &
84!$OMP NUM_THREADS(nthd)
85do ias=1,natmtot
86 is=idxis(ias)
87 nrc=nrcmt(is)
88 nrci=nrcmti(is)
89 irco=nrci+1
90 zfmt(1:npcmt(is))=0.d0
91 do ig=1,ngvc
92 ifg=igfc(ig)
93 do irc=1,nrc
94 x=gc(ig)*rcmt(irc,is)
95 call sbessel(lmax,x,jl(:,irc))
96 end do
97 z1=zfft(ifg)*sfacg(ig,ias)
98 do l=0,lmax
99 do lm=l**2+1,(l+1)**2
100 z2=z1*conjg(ylmg(lm,ig))
101 i1=lmmaxi*(nrci-1)+lm
102 zfmt(lm:i1:lmmaxi)=zfmt(lm:i1:lmmaxi)+jl(l,1:nrci)*z2
103 i0=i1+lmmaxi
104 i1=lmmaxo*(nrc-irco)+i0
105 zfmt(i0:i1:lmmaxo)=zfmt(i0:i1:lmmaxo)+jl(l,irco:nrc)*z2
106 end do
107 end do
108 end do
109 call ztorfmt(nrc,nrci,zfmt,rhomt(:,ias))
110end do
111!$OMP END PARALLEL DO
112call freethd(nthd)
113! convert the density from a coarse to a fine radial mesh
114call rfmtctof(rhomt)
115! add the atomic charge density and the excess charge in each muffin-tin
116t1=chgexs/omega
117do ias=1,natmtot
118 is=idxis(ias)
119 nr=nrmt(is)
120 nri=nrmti(is)
121 iro=nri+1
122 i1=lmmaxi*(nri-1)+1
123 rhomt(1:i1:lmmaxi,ias)=rhomt(1:i1:lmmaxi,ias)+(t1+rhosp(1:nri,is))*y00i
124 i0=i1+lmmaxi
125 i1=lmmaxo*(nr-iro)+i0
126 rhomt(i0:i1:lmmaxo,ias)=rhomt(i0:i1:lmmaxo,ias)+(t1+rhosp(iro:nr,is))*y00i
127end do
128! deallocate the global array rhosp as it is not used again
129deallocate(rhosp)
130! interstitial density determined from the atomic tails and excess charge
131call zfftifc(3,ngdgc,1,zfft)
132do ir=1,ngtc
133 rhoir(ir)=dble(zfft(ir))+t1
134! make sure that the density is always positive
135 if (rhoir(ir) < 1.d-10) rhoir(ir)=1.d-10
136end do
137! convert the interstitial density from coarse to fine grid
138call rfirctof(rhoir,rhoir)
139end subroutine
140!EOC
141
real(8), dimension(:,:), allocatable rhosp
Definition modmain.f90:137
real(8) gmaxvr
Definition modmain.f90:384
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00i
Definition modmain.f90:1234
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
real(8), dimension(:,:), allocatable rcmt
Definition modmain.f90:177
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer lmmaxi
Definition modmain.f90:207
real(8) omega
Definition modmain.f90:20
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
complex(8), dimension(:,:), allocatable sfacg
Definition modmain.f90:430
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) epslat
Definition modmain.f90:24
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8) chgexs
Definition modmain.f90:724
integer, dimension(maxspecies) nrsp
Definition modmain.f90:107
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
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
subroutine rhoinit
Definition rhoinit.f90:10
subroutine sbessel(lmax, x, jl)
Definition sbessel.f90:10
subroutine wsplint(n, x, w)
Definition wsplint.f90:7
subroutine zfftifc(nd, n, sgn, z)
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition ztorfmt.f90:7