The Elk Code
Loading...
Searching...
No Matches
occupyuv.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2019 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
occupyuv
7
use
modmain
8
use
modbog
9
use
modmpi
10
implicit none
11
integer
ik,ist
12
real
(8) chg,w,vn
13
real
(8) e0,e1,e,t1
14
! determine the total charge and fermionic anomalous correlation entropy
15
chg=0.d0
16
face
=0.d0
17
do
ik=1,
nkpt
18
w=
wkpt
(ik)
19
do
ist=1,
nstsv
20
vn=
vnorm
(ist,ik)
21
chg=chg+w*vn
22
if
((vn > 0.d0).and.(vn < 1.d0))
then
23
face
=
face
+w*(vn*log(vn)+(1.d0-vn)*log(1.d0-vn))
24
end if
25
end do
26
end do
27
chg=
occmax
*chg
28
face
=-
occmax
*
face
29
! adjust the Fermi energy
30
efermi
=
efermi
+
tauefm
*(
chgval
-chg)
31
if
(
mp_mpi
)
then
32
if
(abs(chg-
chgval
) >
epschg
)
then
33
write
(*,*)
34
write
(*,
'("Warning(occupyuv): incorrect charge : ",2G18.10)'
) chg,
chgval
35
end if
36
end if
37
! estimate the indirect band gap
38
e0=-1.d8
39
e1=1.d8
40
ikgap
(1)=1
41
ikgap
(2)=1
42
do
ist=1,
nstsv
43
do
ik=1,
nkpt
44
e=
evaluv
(ist,ik)
45
if
(
vnorm
(ist,ik) > 0.5d0) e=-e
46
if
(e <= 0.d0)
then
47
if
(e > e0)
then
48
e0=e
49
ikgap
(1)=ik
50
end if
51
else
52
if
(e < e1)
then
53
e1=e
54
ikgap
(2)=ik
55
end if
56
end if
57
end do
58
end do
59
bandgap
(1)=e1-e0
60
! estimate the direct band gap
61
e=1.d8
62
ikgap
(3)=1
63
do
ik=1,
nkpt
64
e0=-1.d8
65
e1=1.d8
66
do
ist=1,
nstsv
67
t1=
evaluv
(ist,ik)
68
if
(
vnorm
(ist,ik) > 0.5d0) t1=-t1
69
if
(t1 <= 0.d0)
then
70
if
(t1 > e0) e0=t1
71
else
72
if
(t1 < e1) e1=t1
73
end if
74
end do
75
t1=e1-e0
76
if
(t1 < e)
then
77
e=t1
78
ikgap
(3)=ik
79
end if
80
end do
81
bandgap
(2)=e
82
end subroutine
83
modbog
Definition
modbog.f90:6
modbog::tauefm
real(8) tauefm
Definition
modbog.f90:19
modbog::face
real(8) face
Definition
modbog.f90:27
modbog::evaluv
real(8), dimension(:,:), allocatable evaluv
Definition
modbog.f90:15
modbog::vnorm
real(8), dimension(:,:), allocatable vnorm
Definition
modbog.f90:17
modmain
Definition
modmain.f90:6
modmain::wkpt
real(8), dimension(:), allocatable wkpt
Definition
modmain.f90:475
modmain::efermi
real(8) efermi
Definition
modmain.f90:904
modmain::chgval
real(8) chgval
Definition
modmain.f90:722
modmain::ikgap
integer, dimension(3) ikgap
Definition
modmain.f90:914
modmain::bandgap
real(8), dimension(2) bandgap
Definition
modmain.f90:912
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::epschg
real(8) epschg
Definition
modmain.f90:712
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmpi
Definition
modmpi.f90:6
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
occupyuv
subroutine occupyuv
Definition
occupyuv.f90:7
occupyuv.f90
Generated by
1.9.8