The Elk Code
 
Loading...
Searching...
No Matches
rfint.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007 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
6real(8) function rfint(rfmt,rfir)
7use modmain
8implicit none
9! arguments
10real(8), intent(in) :: rfmt(npmtmax,natmtot),rfir(ngtot)
11! local variables
12integer is,ias
13! external functions
14real(8), external :: rfmtint
15! interstitial contribution
16rfint=dot_product(rfir(1:ngtot),cfunir(1:ngtot))
18! muffin-tin contribution
19do ias=1,natmtot
20 is=idxis(ias)
21 rfint=rfint+rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),rfmt(:,ias))
22end do
23end function
24
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer ngtot
Definition modmain.f90:390
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
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
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8) function rfint(rfmt, rfir)
Definition rfint.f90:7
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7