The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine writedosu
7use modmain
8use modulr
9implicit none
10! local variables
11integer nsk(3),ik0,iw
12real(8) dw
13! allocatable arrays
14real(8), allocatable :: w(:),f(:,:),g(:)
15! no k-point reduction
17reducek=0
18! initialise global variables
19call init0
20call init1
21call initulr
22! read in the Fermi energy from STATE_ULR.OUT
23call readstulr
24! get the ULR eigenvalues from file
25do ik0=1,nkpt0
26 call getevalu(ik0)
27end do
28! subtract the Fermi energy
30! generate frequency grid
31allocate(w(nwplot))
32dw=(wplot(2)-wplot(1))/dble(nwplot)
33do iw=1,nwplot
34 w(iw)=dw*dble(iw-1)+wplot(1)
35end do
36! number of subdivisions used for interpolation in the Brillouin zone
37nsk(:)=max(ngrkf/ngridk(:),1)
38! integrate over the Brillouin zone
39allocate(f(nstulr,nkpt0),g(nwplot))
40! normalise DOS to the unit cell
41f(1:nstulr,1:nkpt0)=1.d0/dble(nkpa)
43! write total DOS to file
44open(50,file='TDOSULR.OUT',form='FORMATTED',action='WRITE')
45do iw=1,nwplot
46 write(50,'(2G18.10)') w(iw),g(iw)
47end do
48close(50)
49write(*,*)
50write(*,'("Info(writedosu):")')
51write(*,'(" Ultra long-range total density of states written to TDOSULR.OUT")')
52write(*,*)
53write(*,'(" Fermi energy is at zero in plot")')
54write(*,*)
55write(*,'(" DOS units are states/Hartree/unit cell")')
56deallocate(w,f,g)
57! restore original parameters
59end subroutine
60
subroutine brzint(nsm, ngridk, nsk, ivkik, nw, wint, n, ld, e, f, g)
Definition brzint.f90:10
subroutine getevalu(ik0)
Definition getevalu.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine initulr
Definition initulr.f90:7
real(8) efermi
Definition modmain.f90:904
integer ngrkf
Definition modmain.f90:1072
integer nwplot
Definition modmain.f90:1070
integer nswplot
Definition modmain.f90:1074
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer, dimension(3) ngridk
Definition modmain.f90:448
integer reducek
Definition modmain.f90:455
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer reducek0
Definition modmain.f90:455
integer nstulr
Definition modulr.f90:94
integer nkpa
Definition modulr.f90:24
real(8), dimension(:,:), allocatable evalu
Definition modulr.f90:96
integer nkpt0
Definition modulr.f90:18
subroutine readstulr
Definition readstulr.f90:7
subroutine writedosu
Definition writedosu.f90:7