File: zfinp.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (65 lines) | stat: -rw-r--r-- 2,025 bytes parent folder | download
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
65

! 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: zfinp
! !INTERFACE:
complex(8) function zfinp(zfmt1,zfir1,zfmt2,zfir2)
! !USES:
use modmain
use modomp
! !INPUT/OUTPUT PARAMETERS:
!   zfmt1 : first complex function in spherical harmonics/coordinates for all
!           muffin-tins (in,complex(npcmtmax,natmtot))
!   zfir1 : first complex interstitial function in real-space
!           (in,complex(ngtc))
!   zfmt2 : second complex function in spherical harmonics/coordinates for all
!           muffin-tins (in,complex(npcmtmax,natmtot))
!   zfir2 : second complex interstitial function in real-space
!           (in,complex(ngtc))
! !DESCRIPTION:
!   Calculates the inner product of two complex fuctions over the entire unit
!   cell. The muffin-tin functions should be stored on the coarse radial grid.
!   In the interstitial region, the integrand is multiplied with the
!   characteristic function to remove the contribution from the muffin-tin. See
!   routines {\tt zfmtinp} and {\tt gencfun}.
!
! !REVISION HISTORY:
!   Created July 2004 (Sharma)
!EOP
!BOC
implicit none
! arguments
complex(8), intent(in) :: zfmt1(npcmtmax,natmtot),zfir1(ngtc)
complex(8), intent(in) :: zfmt2(npcmtmax,natmtot),zfir2(ngtc)
! local variables
integer is,ias,ir,nthd
! external functions
complex(8) zfmtinp
external zfmtinp
! interstitial contribution
zfinp=cfrc(1)*conjg(zfir1(1))*zfir2(1)
do ir=2,ngtc
  zfinp=zfinp+cfrc(ir)*conjg(zfir1(ir))*zfir2(ir)
end do
zfinp=zfinp*(omega/dble(ngtc))
! muffin-tin contribution
call omp_hold(natmtot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(is) REDUCTION(+:zfinp) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
  is=idxis(ias)
  zfinp=zfinp+zfmtinp(nrcmt(is),nrcmti(is),rcmt(:,is),r2cmt(:,is), &
   zfmt1(:,ias),zfmt2(:,ias))
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
return
end function
!EOC