The Elk Code
 
Loading...
Searching...
No Matches
chargeu.f90
Go to the documentation of this file.
1
2! Copyright (C) 2019 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
6subroutine chargeu
7use modmain
8use modulr
9use modtest
10implicit none
11! local variables
12integer ifq,is,ias
13integer nrc,nrci
14real(8) t1
15! automatic arrays
16real(8) rfft(nqpt)
17complex(8) zfft(nfqrz)
18! external functions
19real(8), external :: ddot
20complex(8), external :: zfmtint
21! calculate muffin-tin charges
22chgmttot=0.d0
23do ias=1,natmtot
24 is=idxis(ias)
25 nrc=nrcmt(is)
26 nrci=nrcmti(is)
27 do ifq=1,nfqrz
28 zfft(ifq)=zfmtint(nrc,nrci,wr2cmt(:,is),rhoqmt(:,ias,ifq))
29 end do
30 chgmt(ias)=dble(zfft(1))
32 call rzfftifc(3,ngridq,1,rfft,zfft)
33 chgmtru(ias,:)=rfft(:)
34end do
35! calculate interstitial charge
36t1=ddot(ngtc,rhoqir(:,1),2,cfrc,1)
37chgir=t1*omega/dble(ngtc)
38! total calculated charge
40! write muffin-tin charges to file
41call writetest(730,'ULR muffin-tin charges',nv=natmtot*nqpt,tol=5.d-5, &
42 rva=chgmtru)
43end subroutine
44
subroutine chargeu
Definition chargeu.f90:7
real(8) chgmttot
Definition modmain.f90:734
real(8) chgcalc
Definition modmain.f90:728
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
real(8) omega
Definition modmain.f90:20
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer natmtot
Definition modmain.f90:40
real(8), dimension(:), allocatable chgmt
Definition modmain.f90:732
real(8) chgir
Definition modmain.f90:730
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer, dimension(3) ngridq
Definition modmain.f90:515
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
complex(8), dimension(:,:,:), allocatable rhoqmt
Definition modulr.f90:59
real(8), dimension(:,:), allocatable chgmtru
Definition modulr.f90:55
complex(8), dimension(:,:), allocatable rhoqir
Definition modulr.f90:59
subroutine rzfftifc(nd, n, sgn, r, z)
pure complex(8) function zfmtint(nr, nri, wr, zfmt)
Definition zfmtint.f90:7