-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathtetiw.f90
86 lines (59 loc) · 2.13 KB
/
tetiw.f90
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
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
! The code was developed at the Fritz Haber Institute, and
! the intellectual properties and copyright of this file
! are with the Max Planck Society. When you use it, please
! cite R. Gomez-Abal, X. Li, C. Ambrosch-Draxl, M. Scheffler,
! Extended linear tetrahedron method for the calculation of q-dependent
! dynamical response functions, to be published in Comp. Phys. Commun. (2010)
!BOP
!
! !ROUTINE: tetiw
!
! !INTERFACE:
subroutine tetiw(nik,nt,nb,ebd,tetc,wtet,v,efer,iw)
!
! !DESCRIPTION:
!
! This subroutine gives the weight of on one k-point for a certain band
! in an operator integration.
!
! !USES:
use tetra_internal
implicit none
! !INPUT PARAMETERS:
integer(4), intent(in) :: nik ! Number of irreducible
! k-points
integer(4), intent(in) :: nt ! Number of tetrahedra
integer(4), intent(in) :: nb ! Number of bands
real(8), target, intent(in) :: ebd(nb,nik) ! Band energies
integer(4), target, intent(in) :: tetc(4,*) ! id. numbers of
! the corners
! of the tetrahedra
integer(4), target, intent(in) :: wtet(*) ! weight of each
! tetrahedron
real(8), intent(in) :: v ! the volume of the tetrahedra
real(8), intent(in) :: efer ! fermi energy
! !OUTPUT PARAMETERS:
real(8), intent(out) :: iw(nb,nik) ! the value of the integral
! !INTRINSIC ROUTINES:
intrinsic size
! !EXTERNAL ROUTINES:
external intw
! !REVISION HISTORY:
!
! Created: 4th. March 2004 by RGA
!
!EOP
!BOC
nirkp = nik
ntet = nt
nband = nb
vt = v
tetcorn => tetc(1:4,1:ntet)
tetweig => wtet(1:ntet)
eband => ebd(1:nband,1:nik)
!
! intw is called to calculate the weigths
!
call intw(efer,iw)
end subroutine tetiw
!EOC