The Elk Code
 
Loading...
Searching...
No Matches
charge.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: charge
8! !INTERFACE:
9subroutine charge
10! !USES:
11use modmain
12use modtest
13! !DESCRIPTION:
14! Computes the muffin-tin, interstitial and total charges by integrating the
15! density.
16!
17! !REVISION HISTORY:
18! Created April 2003 (JKD)
19!EOP
20!BOC
21implicit none
22! local variables
23integer is,ias
24! external functions
25real(8), external :: rfmtint
26! find the muffin-tin charges
27chgmttot=0.d0
28do ias=1,natmtot
29 is=idxis(ias)
30 chgmt(ias)=rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),rhomt(:,ias))
32end do
33! find the interstitial charge
34chgir=(omega/ngtot)*dot_product(rhoir(1:ngtot),cfunir(1:ngtot))
35! total calculated charge
37! write total calculated charge to test file
38call writetest(400,'calculated total charge',tol=1.d-6,rv=chgcalc)
39end subroutine
40!EOC
41
subroutine charge
Definition charge.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8) chgmttot
Definition modmain.f90:734
integer ngtot
Definition modmain.f90:390
real(8) chgcalc
Definition modmain.f90:728
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
real(8) omega
Definition modmain.f90:20
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
real(8), dimension(:), allocatable cfunir
Definition modmain.f90:436
real(8), dimension(:), allocatable chgmt
Definition modmain.f90:732
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8) chgir
Definition modmain.f90:730
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7