-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathdostet.f90
76 lines (54 loc) · 1.78 KB
/
dostet.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
! 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: dostet
!
! !INTERFACE:
real(8) function dostet(nik,nbd,eband,ntet,tetc,wtet,vt,e)
!
! !DESCRIPTION:
!
! This subroutine calculates the density of states at an energy e
!
!
! !USES:
use order
implicit none
! !INPUT PARAMETERS:
integer(4), intent(in) :: nik ! Number of irreducible k-points
integer(4), intent(in) :: nbd ! Maximum number of bands
real(8), intent(in) :: eband(nbd,nik) ! Band energies
integer(4), intent(in) :: ntet ! Number of tetrahedra
integer(4), intent(in) :: tetc(4,*)! id. numbers of the corners
! of the tetrahedra
integer(4), intent(in) :: wtet(*) ! weight of each tetrahedron
real(8), intent(in) :: vt ! the volume of the tetrahedra
real(8), intent(in) :: e ! energy
! !LOCAL VARIABLES:
integer(4) :: itet,i,ib
real(8), dimension(4) :: ee
real(8), external :: dos1t
! !REVISION HISTORY:
!
! Created: 3th. March 2004. by RGA
!
!EOP
!
!BOC
dostet=0.0d0
do itet=1,ntet
do ib=1,nbd
do i=1,4
ee(i)=eband(ib,tetc(i,itet))
enddo
call sort(4,ee)
dostet=dostet+wtet(itet)*dos1t(ee,e,vt)
enddo
enddo
end function dostet
!EOC