The Elk Code
writedosu.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2024 Wenhan Chen, J. K. Dewhurst and S. Sharma.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
writedosu
7
use
modmain
8
use
modulr
9
implicit none
10
! local variables
11
integer
nsk(3),ik0,iw
12
real(8)
dw
13
! allocatable arrays
14
real(8)
,
allocatable
:: w(:),f(:,:),g(:)
15
! no k-point reduction
16
reducek0
=
reducek
17
reducek
=0
18
! initialise global variables
19
call
init0
20
call
init1
21
call
initulr
22
! read in the Fermi energy from STATE_ULR.OUT
23
call
readstulr
24
! get the ULR eigenvalues from file
25
do
ik0=1,
nkpt0
26
call
getevalu
(ik0)
27
end do
28
! subtract the Fermi energy
29
evalu
(1:
nstulr
,1:
nkpt0
)=
evalu
(1:
nstulr
,1:
nkpt0
)-
efermi
30
! generate frequency grid
31
allocate
(w(
nwplot
))
32
dw=(
wplot
(2)-
wplot
(1))/dble(
nwplot
)
33
do
iw=1,
nwplot
34
w(iw)=dw*dble(iw-1)+
wplot
(1)
35
end do
36
! number of subdivisions used for interpolation in the Brillouin zone
37
nsk(:)=max(
ngrkf
/
ngridk
(:),1)
38
! integrate over the Brillouin zone
39
allocate
(f(
nstulr
,
nkpt0
),g(
nwplot
))
40
! normalise DOS to the unit cell
41
f(1:
nstulr
,1:
nkpt0
)=1.d0/dble(
nkpa
)
42
call
brzint
(
nswplot
,
ngridk
,nsk,
ivkik
,
nwplot
,
wplot
,
nstulr
,
nstulr
,
evalu
,f,g)
43
! write total DOS to file
44
open
(50,file=
'TDOSULR.OUT'
,form=
'FORMATTED'
,action=
'WRITE'
)
45
do
iw=1,
nwplot
46
write
(50,
'(2G18.10)'
) w(iw),g(iw)
47
end do
48
close
(50)
49
write
(*,*)
50
write
(*,
'("Info(writedosu):")'
)
51
write
(*,
'(" Ultra long-range total density of states written to TDOSULR.OUT")'
)
52
write
(*,*)
53
write
(*,
'(" Fermi energy is at zero in plot")'
)
54
write
(*,*)
55
write
(*,
'(" DOS units are states/Hartree/unit cell")'
)
56
deallocate
(w,f,g)
57
! restore original parameters
58
reducek
=
reducek0
59
end subroutine
60
readstulr
subroutine readstulr
Definition:
readstulr.f90:7
modmain::ngrkf
integer ngrkf
Definition:
modmain.f90:1075
modmain::efermi
real(8) efermi
Definition:
modmain.f90:907
modulr::nstulr
integer nstulr
Definition:
modulr.f90:95
modmain::reducek0
integer reducek0
Definition:
modmain.f90:455
modulr::nkpt0
integer nkpt0
Definition:
modulr.f90:18
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition:
modmain.f90:467
brzint
subroutine brzint(nsm, ngridk, nsk, ivkik, nw, wint, n, ld, e, f, g)
Definition:
brzint.f90:10
modmain
Definition:
modmain.f90:6
modmain::wplot
real(8), dimension(2) wplot
Definition:
modmain.f90:1079
modmain::ngridk
integer, dimension(3) ngridk
Definition:
modmain.f90:448
init1
subroutine init1
Definition:
init1.f90:10
writedosu
subroutine writedosu
Definition:
writedosu.f90:7
modmain::nswplot
integer nswplot
Definition:
modmain.f90:1077
modmain::reducek
integer reducek
Definition:
modmain.f90:455
init0
subroutine init0
Definition:
init0.f90:10
getevalu
subroutine getevalu(ik0)
Definition:
getevalu.f90:7
modulr
Definition:
modulr.f90:6
modmain::nwplot
integer nwplot
Definition:
modmain.f90:1073
modulr::nkpa
integer nkpa
Definition:
modulr.f90:24
modulr::evalu
real(8), dimension(:,:), allocatable evalu
Definition:
modulr.f90:97
initulr
subroutine initulr
Definition:
initulr.f90:7
writedosu.f90
Generated by
1.8.14