File: rhonorm.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (71 lines) | stat: -rw-r--r-- 1,796 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: rhonorm
! !INTERFACE:
subroutine rhonorm
! !USES:
use modmain
! !DESCRIPTION:
!   Loss of precision of the calculated total charge can result because the
!   muffin-tin density is computed on a set of $(\theta,\phi)$ points and then
!   transformed to a spherical harmonic representation. This routine adds a
!   constant to the density so that the total charge is correct. If the error in
!   total charge exceeds a certain tolerance then a warning is issued.
!
! !REVISION HISTORY:
!   Created April 2003 (JKD)
!   Changed from rescaling to adding, September 2006 (JKD)
!EOP
!BOC
implicit none
! local variables
integer is,ia,ias
integer nr,nri,ir,i
real(8) t1,t2
if (.not.trhonorm) return
! check error in total charge
t1=chgcalc/chgtot-1.d0
if (abs(t1).gt.epschg) then
  write(*,*)
  write(*,'("Warning(rhonorm): total charge density incorrect for s.c. &
   &loop ",I5)') iscl
  write(*,'(" Calculated : ",G18.10)') chgcalc
  write(*,'(" Required   : ",G18.10)') chgtot
end if
! error in average density
t1=(chgtot-chgcalc)/omega
! add the constant difference to the density
t2=t1/y00
do ias=1,natmtot
  is=idxis(ias)
  nr=nrmt(is)
  nri=nrmti(is)
  i=1
  do ir=1,nri
    rhomt(i,ias)=rhomt(i,ias)+t2
    i=i+lmmaxi
  end do
  do ir=nri+1,nr
    rhomt(i,ias)=rhomt(i,ias)+t2
    i=i+lmmaxo
  end do
end do
rhoir(:)=rhoir(:)+t1
! add the difference to the charges
do is=1,nspecies
  t2=t1*(4.d0*pi/3.d0)*rmt(is)**3
  do ia=1,natoms(is)
    ias=idxas(ia,is)
    chgmt(ias)=chgmt(ias)+t2
    chgmttot=chgmttot+t2
  end do
end do
chgir=chgtot-chgmttot
return
end subroutine
!EOC