The Elk Code
 
Loading...
Searching...
No Matches
occupy.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: occupy
8! !INTERFACE:
9subroutine occupy
10! !USES:
11use modmain
12use modtest
13! !DESCRIPTION:
14! Finds the Fermi energy and sets the occupation numbers for the
15! second-variational states using the routine {\tt fermi}.
16!
17! !REVISION HISTORY:
18! Created February 2004 (JKD)
19! Added gap estimation, November 2009 (F. Cricchio)
20! Added adaptive smearing width, April 2010 (T. Bjorkman)
21!EOP
22!BOC
23implicit none
24! local variables
25integer, parameter :: maxit=1000
26integer ik,ist,it
27real(8) e0,e1,e,de
28real(8) ei0,ei1,ed0,ed1
29real(8) chg,w,x,t0
30! external functions
31real(8), external :: sdelta,stheta
32! determine the smearing width automatically if required
33if ((autoswidth).and.(iscl > 1)) call findswidth
34! find minimum and maximum eigenvalues
35e0=evalsv(1,1)
36e1=e0
37do ik=1,nkpt
38 do ist=1,nstsv
39 e=evalsv(ist,ik)
40 if (e < e0) e0=e
41 if (e > e1) e1=e
42 end do
43end do
44if (e0 < e0min) then
45 write(*,*)
46 write(*,'("Warning(occupy): minimum eigenvalue less than minimum &
47 &linearisation energy : ",2G18.10)') e0,e0min
48 write(*,'(" for s.c. loop ",I5)') iscl
49end if
50t0=1.d0/swidth
51! determine the Fermi energy using the bisection method
52do it=1,maxit
53 efermi=0.5d0*(e0+e1)
54 chg=0.d0
55 do ik=1,nkpt
56 w=wkpt(ik)
57 do ist=1,nstsv
58 e=evalsv(ist,ik)
59 if (e < e0min) then
60 occsv(ist,ik)=0.d0
61 else
62 x=(efermi-e)*t0
63 occsv(ist,ik)=occmax*stheta(stype,x)
64 chg=chg+w*occsv(ist,ik)
65 end if
66 end do
67 end do
68 if (chg < chgval) then
69 e0=efermi
70 else
71 e1=efermi
72 end if
73 if ((e1-e0) < 1.d-12) goto 10
74end do
75write(*,*)
76write(*,'("Warning(occupy): could not find Fermi energy")')
7710 continue
78! find the density of states at the Fermi surface in units of
79! states/Hartree/unit cell
80fermidos=0.d0
81do ik=1,nkpt
82 w=wkpt(ik)
83 do ist=1,nstsv
84 x=(evalsv(ist,ik)-efermi)*t0
86 end do
87 if (abs(occsv(nstsv,ik)) > epsocc) then
88 write(*,*)
89 write(*,'("Warning(occupy): not enough empty states for k-point ",I6)') ik
90 write(*,'(" and s.c. loop ",I5)') iscl
91 end if
92end do
94! write Fermi density of states to test file
95call writetest(500,'DOS at Fermi energy',tol=5.d-3,rv=fermidos)
96! estimate the indirect and direct band gaps (FC)
97ei0=-1.d8; ei1=1.d8
98de=1.d8
99ikgap(1:3)=1
100do ik=1,nkpt
101 ed0=-1.d8; ed1=1.d8
102 do ist=1,nstsv
103 e=evalsv(ist,ik)
104 if (e <= efermi) then
105 if (e > ed0) ed0=e
106 if (e > ei0) then
107! transfer is a workaround for a bug in Intel Fortran versions 17 and 18
108 ikgap(1)=transfer(ik,ik)
109 ei0=e
110 end if
111 else
112 if (e < ed1) ed1=e
113 if (e < ei1) then
114 ikgap(2)=ik
115 ei1=e
116 end if
117 end if
118 end do
119 e=ed1-ed0
120 if (e < de) then
121 ikgap(3)=ik
122 de=e
123 end if
124end do
125bandgap(1)=ei1-ei0
126bandgap(2)=de
127! write band gap to test file
128call writetest(510,'estimated indirect band gap',tol=2.d-2,rv=bandgap(1))
129end subroutine
130!EOC
131
subroutine findswidth
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
logical autoswidth
Definition modmain.f90:894
real(8) efermi
Definition modmain.f90:904
real(8) chgval
Definition modmain.f90:722
integer, dimension(3) ikgap
Definition modmain.f90:914
real(8), dimension(2) bandgap
Definition modmain.f90:912
real(8) e0min
Definition modmain.f90:825
real(8) swidth
Definition modmain.f90:892
integer nkpt
Definition modmain.f90:461
real(8) fermidos
Definition modmain.f90:910
integer nstsv
Definition modmain.f90:886
real(8) epsocc
Definition modmain.f90:900
real(8) occmax
Definition modmain.f90:898
integer stype
Definition modmain.f90:888
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer iscl
Definition modmain.f90:1048
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine occupy
Definition occupy.f90:10
real(8) function sdelta(stype, x)
Definition sdelta.f90:10
real(8) function stheta(stype, x)
Definition stheta.f90:10