The Elk Code
 
Loading...
Searching...
No Matches
rfinp.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: rfinp
8! !INTERFACE:
9real(8) function rfinp(rfmt1,rfir1,rfmt2,rfir2)
10! !USES:
11use modmain
12use modomp
13! !INPUT/OUTPUT PARAMETERS:
14! rfmt1 : first function in real spherical harmonics for all muffin-tins
15! (in,real(npmtmax,natmtot))
16! rfir1 : first real interstitial function in real-space (in,real(ngtot))
17! rfmt2 : second function in real spherical harmonics for all muffin-tins
18! (in,real(npmtmax,natmtot))
19! rfir2 : second real interstitial function in real-space (in,real(ngtot))
20! !DESCRIPTION:
21! Calculates the inner product of two real functions over the entire unit
22! cell. The input muffin-tin functions should have angular momentum cut-off
23! {\tt lmaxo}. In the interstitial region, the integrand is multiplied with
24! the characteristic function, $\tilde{\Theta}({\bf r})$, to remove the
25! contribution from the muffin-tin. See routines {\tt rfmtinp} and
26! {\tt gencfun}.
27!
28! !REVISION HISTORY:
29! Created July 2004 (JKD)
30!EOP
31!BOC
32implicit none
33! arguments
34real(8), intent(in) :: rfmt1(npmtmax,natmtot),rfir1(ngtot)
35real(8), intent(in) :: rfmt2(npmtmax,natmtot),rfir2(ngtot)
36! local variables
37integer is,ias,nthd
38real(8) t1
39! external functions
40real(8), external :: rfmtinp
41rfinp=0.d0
42call holdthd(natmtot+1,nthd)
43!$OMP PARALLEL DEFAULT(SHARED) &
44!$OMP PRIVATE(is,t1) REDUCTION(+:rfinp) &
45!$OMP NUM_THREADS(nthd)
46!$OMP DO SCHEDULE(DYNAMIC)
47do ias=1,natmtot
48 is=idxis(ias)
49! muffin-tin contribution
50 rfinp=rfinp+rfmtinp(nrmt(is),nrmti(is),wr2mt(:,is),rfmt1(:,ias),rfmt2(:,ias))
51end do
52!$OMP END DO NOWAIT
53! interstitial contribution
54!$OMP SINGLE
55t1=sum(rfir1(1:ngtot)*rfir2(1:ngtot)*cfunir(1:ngtot))
57!$OMP END SINGLE
58!$OMP END PARALLEL
59call freethd(nthd)
60end function
61!EOC
62
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
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
real(8) function rfinp(rfmt1, rfir1, rfmt2, rfir2)
Definition rfinp.f90:10
pure real(8) function rfmtinp(nr, nri, wr, rfmt1, rfmt2)
Definition rfmtinp.f90:10