The Elk Code
 
Loading...
Searching...
No Matches
ephdos.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3! E. K. U. Gross. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine ephdos
7use modmain
8use modphonon
9use modbog
10implicit none
11! local variables
12integer nsk(3),ik,ist,iw
13real(8) v(3),dw,vn,t1
14! allocatable arrays
15integer, allocatable :: idx(:)
16real(8), allocatable :: w(:),f(:,:),g(:)
17! initialise universal variables
18call init0
19call init1
20call init2
21call initeph
22allocate(idx(nstsv),w(nstsv))
23do ik=1,nkpt
24! get the eigenvalues from file
25 call getevaluv(ik,evaluv(:,ik))
26! negate eigenvalues corresponding to V-norm > 1/2
27 do ist=1,nstsv
28 if (vnorm(ist,ik) > 0.5d0) evaluv(ist,ik)=-evaluv(ist,ik)
29 end do
30! arrange in ascending order
31 w(:)=evaluv(:,ik)
32 call sortidx(nstsv,w,idx)
33 evaluv(:,ik)=w(idx(:))
34! put the V-norm into the same order as the eigenvalues
35 w(:)=vnorm(:,ik)
36 vnorm(:,ik)=w(idx(:))
37end do
38deallocate(idx,w)
39! generate the partial and total DOS and write to file
40allocate(w(nwplot),f(nstsv,nkptnr),g(nwplot))
41! generate frequency grid
42dw=(wplot(2)-wplot(1))/dble(nwplot)
43do iw=1,nwplot
44 w(iw)=dw*dble(iw-1)+wplot(1)
45end do
46! number of subdivisions used for interpolation in the Brillouin zone
47nsk(:)=max(ngrkf/ngridk(:),1)
48! set the weight array
49f(:,:)=occmax
50! integrate over the Brillouin zone
52! output the total electronic DOS to file
53open(50,file='TDOS_EPH.OUT.OUT',form='FORMATTED',action='WRITE')
54do iw=1,nwplot
55 write(50,'(2G18.10)') w(iw),g(iw)
56end do
57close(50)
58! output the FACE vs energy histogram to file
59open(50,file='FACEEH.OUT',form='FORMATTED',action='WRITE')
60do ik=1,nkpt
61! map k-vector to first Brillouin zone
62 v(:)=vkc(:,ik)
63 call vecfbz(epslat,bvec,v)
64 do ist=1,nstsv
65 vn=vnorm(ist,ik)
66 if ((vn > 0.d0).and.(vn < 1.d0)) then
67 t1=-(vn*log(vn)+(1.d0-vn)*log(1.d0-vn))
68 else
69 t1=0.d0
70 end if
71 if (t1 < 1.d-4) cycle
72 write(50,'(5G18.10)') evaluv(ist,ik),t1,v
73 end do
74end do
75close(50)
76write(*,*)
77write(*,'("Info(ephdos):")')
78write(*,'(" Total electronic density of states for the electron-phonon")')
79write(*,'(" system written to TDOS_EPH.OUT.OUT")')
80write(*,*)
81write(*,'(" Fermionic anomalous correlation entropy vs energy histogram")')
82write(*,'(" written to FACEEH.OUT")')
83write(*,*)
84write(*,'(" Fermi energy is at zero in plots")')
85write(*,*)
86write(*,'(" DOS units are states/Hartree/unit cell")')
87deallocate(w,f,g)
88end subroutine
89
subroutine brzint(nsm, ngridk, nsk, ivkik, nw, wint, n, ld, e, f, g)
Definition brzint.f90:10
subroutine ephdos
Definition ephdos.f90:7
subroutine getevaluv(ik, evaluvp)
Definition getevaluv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine initeph
Definition initeph.f90:7
real(8), dimension(:,:), allocatable evaluv
Definition modbog.f90:15
real(8), dimension(:,:), allocatable vnorm
Definition modbog.f90:17
integer nkptnr
Definition modmain.f90:463
integer ngrkf
Definition modmain.f90:1072
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer nwplot
Definition modmain.f90:1070
integer nswplot
Definition modmain.f90:1074
real(8) epslat
Definition modmain.f90:24
integer nkpt
Definition modmain.f90:461
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer, dimension(3) ngridk
Definition modmain.f90:448
integer nstsv
Definition modmain.f90:886
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
pure subroutine sortidx(n, x, idx)
Definition sortidx.f90:10
subroutine vecfbz(eps, bvec, vpl)
Definition vecfbz.f90:10