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
|
! 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: rfinp
! !INTERFACE:
real(8) function rfinp(rfmt1,rfir1,rfmt2,rfir2)
! !USES:
use modmain
use modomp
! !INPUT/OUTPUT PARAMETERS:
! rfmt1 : first function in real spherical harmonics for all muffin-tins
! (in,real(npmtmax,natmtot))
! rfir1 : first real interstitial function in real-space (in,real(ngtot))
! rfmt2 : second function in real spherical harmonics for all muffin-tins
! (in,real(npmtmax,natmtot))
! rfir2 : second real interstitial function in real-space (in,real(ngtot))
! !DESCRIPTION:
! Calculates the inner product of two real functions over the entire unit
! cell. The input muffin-tin functions should have angular momentum cut-off
! {\tt lmaxo}. In the interstitial region, the integrand is multiplied with
! the characteristic function, $\tilde{\Theta}({\bf r})$, to remove the
! contribution from the muffin-tin. See routines {\tt rfmtinp} and
! {\tt gencfun}.
!
! !REVISION HISTORY:
! Created July 2004 (JKD)
!EOP
!BOC
implicit none
! arguments
real(8), intent(in) :: rfmt1(npmtmax,natmtot),rfir1(ngtot)
real(8), intent(in) :: rfmt2(npmtmax,natmtot),rfir2(ngtot)
! local variables
integer is,ias,ir,nthd
! external functions
real(8) rfmtinp
external rfmtinp
! interstitial contribution
rfinp=rfir1(1)*rfir2(1)*cfunir(1)
do ir=2,ngtot
rfinp=rfinp+rfir1(ir)*rfir2(ir)*cfunir(ir)
end do
rfinp=rfinp*omega/dble(ngtot)
! muffin-tin contribution
call omp_hold(natmtot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(is) REDUCTION(+:rfinp) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
is=idxis(ias)
rfinp=rfinp+rfmtinp(nrmt(is),nrmti(is),rsp(:,is),r2sp(:,is),rfmt1(:,ias), &
rfmt2(:,ias))
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
return
end function
!EOC
|