The Elk Code
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 ! automatic arrays
33 integer idx(nstsv*nkpt)
34 n=nstsv*nkpt
35 ! find the index which sorts the eigenvalues in ascending order
36 call sortidx(n,evalsv,idx)
37 ! find the highest eigenvalue < efermi
38 m=n
39 e=efermi
40 do i=n,1,-1
41  j=idx(i)
42  ik=(j-1)/nstsv+1
43  ist=mod(j-1,nstsv)+1
44  if (evalsv(ist,ik) < efermi) then
45  m=i
46  e=evalsv(ist,ik)
47  exit
48  end if
49 end do
50 ! search downwards until a gap larger than eps is found
51 do i=m,1,-1
52  j=idx(i)
53  ik=(j-1)/nstsv+1
54  ist=mod(j-1,nstsv)+1
55  if ((e-evalsv(ist,ik)) > eps) exit
56  e=evalsv(ist,ik)
57 end do
58 ! valence band width
59 vbw=efermi-e
60 vbw=max(vbw,1.d-2)
61 ! determine swidth
62 swidth=(sqrt(2.d0*vbw)/mstar)*(6.d0*pi**2/(omega*dble(nkptnr)))**(1.d0/3.d0)
63 end subroutine
64 !EOC
65 
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
real(8) omega
Definition: modmain.f90:20
subroutine findswidth
Definition: findswidth.f90:10
real(8) swidth
Definition: modmain.f90:895
integer nkptnr
Definition: modmain.f90:463
real(8), parameter pi
Definition: modmain.f90:1232
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
real(8) mstar
Definition: modmain.f90:899