The Elk Code
 
Loading...
Searching...
No Matches
rhonorm.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: rhonorm
8! !INTERFACE:
9subroutine rhonorm
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Loss of precision of the calculated total charge can result because the
14! muffin-tin density is computed on a set of $(\theta,\phi)$ points and then
15! transformed to a spherical harmonic representation. This routine adds a
16! constant to the density so that the total charge is correct. If the error in
17! total charge exceeds a certain tolerance then a warning is issued.
18!
19! !REVISION HISTORY:
20! Created April 2003 (JKD)
21! Changed from rescaling to adding, September 2006 (JKD)
22!EOP
23!BOC
24implicit none
25! local variables
26integer is,ia,ias
27integer nr,nri,iro,i0,i1
28real(8) t1,t2
29if (.not.trhonorm) return
30! check error in total charge
31t1=chgcalc/chgtot-1.d0
32if (abs(t1) > epschg) then
33 write(*,*)
34 write(*,'("Warning(rhonorm): total charge density incorrect for s.c. &
35 &loop ",I5)') iscl
36 write(*,'(" Calculated : ",G18.10)') chgcalc
37 write(*,'(" Required : ",G18.10)') chgtot
38end if
39! error in average density
41! add the constant difference to the density
42t2=t1*y00i
43do ias=1,natmtot
44 is=idxis(ias)
45 nr=nrmt(is)
46 nri=nrmti(is)
47 iro=nri+1
48 i1=lmmaxi*(nri-1)+1
49 rhomt(1:i1:lmmaxi,ias)=rhomt(1:i1:lmmaxi,ias)+t2
50 i0=i1+lmmaxi
51 i1=lmmaxo*(nr-iro)+i0
52 rhomt(i0:i1:lmmaxo,ias)=rhomt(i0:i1:lmmaxo,ias)+t2
53end do
54rhoir(1:ngtot)=rhoir(1:ngtot)+t1
55! add the difference to the charges
56t1=t1*(fourpi/3.d0)
57do is=1,nspecies
58 t2=t1*rmt(is)**3
59 do ia=1,natoms(is)
60 ias=idxas(ia,is)
61 chgmt(ias)=chgmt(ias)+t2
63 end do
64end do
66end subroutine
67!EOC
68
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8) chgmttot
Definition modmain.f90:734
real(8), parameter y00i
Definition modmain.f90:1234
integer ngtot
Definition modmain.f90:390
real(8) chgcalc
Definition modmain.f90:728
logical trhonorm
Definition modmain.f90:618
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(maxspecies) rmt
Definition modmain.f90:162
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer lmmaxi
Definition modmain.f90:207
real(8) omega
Definition modmain.f90:20
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) chgtot
Definition modmain.f90:726
integer natmtot
Definition modmain.f90:40
real(8), dimension(:), allocatable chgmt
Definition modmain.f90:732
real(8) epschg
Definition modmain.f90:712
integer lmmaxo
Definition modmain.f90:203
real(8) chgir
Definition modmain.f90:730
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
integer nspecies
Definition modmain.f90:34
integer iscl
Definition modmain.f90:1048
subroutine rhonorm
Definition rhonorm.f90:10