The Elk Code
Loading...
Searching...
No Matches
findswidth.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2010 T. Bjorkman and O. Granas.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: findswidth
8
! !INTERFACE:
9
subroutine
findswidth
10
! !USES:
11
use
modmain
12
! !DESCRIPTION:
13
! Calculates the smearing width from the $k$-point density,
14
! $V_{\rm BZ}/n_k$; the valence band width, $W$; and an effective mass
15
! parameter, $m^{*}$; according to
16
! $$ \sigma=\frac{\sqrt{2W}}{m^{*}}\left(\frac{3}{4\pi}
17
! \frac{V_{\rm BZ}}{n_k}\right)^{1/3}. $$
18
! The valence bandwidth is determined by stepping down in energy from the
19
! Fermi level until a gap larger than a given tolerance is found. This method
20
! was presented in T. Bj\"{o}rkman and O. Gr\aa n\"{a}s, {\it Int. J. Quant.
21
! Chem.} DOI: 10.1002/qua.22476.
22
!
23
! !REVISION HISTORY:
24
! Created April 2010 (Torbjorn Bjorkman and JKD)
25
!EOP
26
!BOC
27
implicit none
28
! local variables
29
integer
i,j,m,n,ik,ist
30
real
(8),
parameter
:: eps=0.05d0
31
real
(8) e,vbw
32
! allocatable arrays
33
integer
,
allocatable
:: idx(:)
34
n=
nstsv
*
nkpt
35
allocate
(idx(n))
36
! find the index which sorts the eigenvalues in ascending order
37
call
sortidx
(n,
evalsv
,idx)
38
! find the highest eigenvalue < efermi
39
m=n
40
e=
efermi
41
do
i=n,1,-1
42
j=idx(i)
43
ik=(j-1)/
nstsv
+1
44
ist=mod(j-1,
nstsv
)+1
45
if
(
evalsv
(ist,ik) <
efermi
)
then
46
m=i
47
e=
evalsv
(ist,ik)
48
exit
49
end if
50
end do
51
! search downwards until a gap larger than eps is found
52
do
i=m,1,-1
53
j=idx(i)
54
ik=(j-1)/
nstsv
+1
55
ist=mod(j-1,
nstsv
)+1
56
if
((e-
evalsv
(ist,ik)) > eps)
exit
57
e=
evalsv
(ist,ik)
58
end do
59
! valence band width
60
vbw=
efermi
-e
61
vbw=max(vbw,1.d-2)
62
! determine swidth
63
swidth
=(sqrt(2.d0*vbw)/
mstar
)*(6.d0*
pi
**2/(
omega
*dble(
nkptnr
)))**(1.d0/3.d0)
64
deallocate
(idx)
65
end subroutine
66
!EOC
67
findswidth
subroutine findswidth
Definition
findswidth.f90:10
modmain
Definition
modmain.f90:6
modmain::efermi
real(8) efermi
Definition
modmain.f90:904
modmain::nkptnr
integer nkptnr
Definition
modmain.f90:463
modmain::pi
real(8), parameter pi
Definition
modmain.f90:1229
modmain::omega
real(8) omega
Definition
modmain.f90:20
modmain::mstar
real(8) mstar
Definition
modmain.f90:896
modmain::swidth
real(8) swidth
Definition
modmain.f90:892
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
sortidx
pure subroutine sortidx(n, x, idx)
Definition
sortidx.f90:10
findswidth.f90
Generated by
1.9.8