The Elk Code
Loading...
Searching...
No Matches
occupyulr.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2018 T. Mueller, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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
occupyulr
7
use
modmain
8
use
modulr
9
implicit none
10
! local variables
11
integer
,
parameter
:: maxit=1000
12
integer
ik0,ik,ist,it
13
real
(8) e0,e1,e
14
real
(8) chg,x,t1
15
! external functions
16
real
(8),
external
:: stheta
17
! find minimum and maximum eigenvalues
18
e0=
evalu
(1,1)
19
e1=e0
20
do
ik0=1,
nkpt0
21
do
ist=1,
nstulr
22
e=
evalu
(ist,ik0)
23
if
(e < e0) e0=e
24
if
(e > e1) e1=e
25
end do
26
end do
27
if
(e0 <
e0min
)
then
28
write
(*,*)
29
write
(*,
'("Warning(occupyulr): minimum eigenvalue less than minimum &
30
&linearisation energy : ",2G18.10)'
) e0,
e0min
31
write
(*,
'(" for s.c. loop ",I5)'
)
iscl
32
end if
33
t1=1.d0/
swidth
34
! determine the Fermi energy using the bisection method
35
do
it=1,maxit
36
efermi
=0.5d0*(e0+e1)
37
chg=0.d0
38
do
ik0=1,
nkpt0
39
! central k-point
40
ik=(ik0-1)*
nkpa
+1
41
do
ist=1,
nstulr
42
e=
evalu
(ist,ik0)
43
if
(e <
e0min
)
then
44
occulr
(ist,ik0)=0.d0
45
else
46
x=(
efermi
-e)*t1
47
occulr
(ist,ik0)=
occmax
*
stheta
(
stype
,x)
48
chg=chg+
wkpt
(ik)*
occulr
(ist,ik0)
49
end if
50
end do
51
end do
52
if
(chg <
chgval
)
then
53
e0=
efermi
54
else
55
e1=
efermi
56
end if
57
if
((e1-e0) < 1.d-12)
return
58
end do
59
write
(*,*)
60
write
(*,
'("Warning(occupyulr): could not find Fermi energy")'
)
61
end subroutine
62
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::e0min
real(8) e0min
Definition
modmain.f90:825
modmain::swidth
real(8) swidth
Definition
modmain.f90:892
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::stype
integer stype
Definition
modmain.f90:888
modmain::iscl
integer iscl
Definition
modmain.f90:1048
modulr
Definition
modulr.f90:6
modulr::nstulr
integer nstulr
Definition
modulr.f90:94
modulr::nkpa
integer nkpa
Definition
modulr.f90:24
modulr::evalu
real(8), dimension(:,:), allocatable evalu
Definition
modulr.f90:96
modulr::occulr
real(8), dimension(:,:), allocatable occulr
Definition
modulr.f90:98
modulr::nkpt0
integer nkpt0
Definition
modulr.f90:18
occupyulr
subroutine occupyulr
Definition
occupyulr.f90:7
stheta
real(8) function stheta(stype, x)
Definition
stheta.f90:10
occupyulr.f90
Generated by
1.9.8