The Elk Code
Loading...
Searching...
No Matches
doccupy.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2015 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
doccupy
7
use
modmain
8
use
modphonon
9
implicit none
10
! local variables
11
integer
,
parameter
:: maxit=1000
12
integer
ik,jk,ist,it
13
real
(8) de0,de1,de
14
real
(8) dchg,x,dx,t1
15
! external functions
16
real
(8),
external
:: sdelta
17
if
(.not.
tphq0
)
return
18
de0=1.d6
19
de1=-1.d6
20
do
ik=1,
nkptnr
21
do
ist=1,
nstsv
22
de=
devalsv
(ist,ik)
23
if
(de < de0) de0=de
24
if
(de > de1) de1=de
25
end do
26
end do
27
t1=1.d0/
swidth
28
do
it=1,maxit
29
defermi
=0.5d0*(de0+de1)
30
dchg=0.d0
31
do
ik=1,
nkptnr
32
jk=
ivkik
(
ivk
(1,ik),
ivk
(2,ik),
ivk
(3,ik))
33
do
ist=1,
nstsv
34
x=(
efermi
-
evalsv
(ist,jk))*t1
35
dx=(
defermi
-
devalsv
(ist,ik))*t1
36
doccsv
(ist,ik)=
occmax
*
sdelta
(
stype
,x)*dx
37
dchg=dchg+
wkptnr
*
doccsv
(ist,ik)
38
end do
39
end do
40
if
(dchg < 0.d0)
then
41
de0=
defermi
42
else
43
de1=
defermi
44
end if
45
if
((de1-de0) < 1.d-12)
return
46
end do
47
write
(*,*)
48
write
(*,
'("Warning(doccupy): could not find Fermi energy derivative")'
)
49
write
(*,
'(" for s.c. loop ",I5)'
)
iscl
50
end subroutine
51
doccupy
subroutine doccupy
Definition
doccupy.f90:7
modmain
Definition
modmain.f90:6
modmain::efermi
real(8) efermi
Definition
modmain.f90:904
modmain::nkptnr
integer nkptnr
Definition
modmain.f90:463
modmain::ivk
integer, dimension(:,:), allocatable ivk
Definition
modmain.f90:465
modmain::swidth
real(8) swidth
Definition
modmain.f90:892
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition
modmain.f90:467
modmain::stype
integer stype
Definition
modmain.f90:888
modmain::iscl
integer iscl
Definition
modmain.f90:1048
modmain::wkptnr
real(8) wkptnr
Definition
modmain.f90:477
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
modphonon
Definition
modphonon.f90:6
modphonon::doccsv
real(8), dimension(:,:), allocatable doccsv
Definition
modphonon.f90:134
modphonon::tphq0
logical tphq0
Definition
modphonon.f90:17
modphonon::devalsv
real(8), dimension(:,:), allocatable devalsv
Definition
modphonon.f90:132
modphonon::defermi
real(8) defermi
Definition
modphonon.f90:128
sdelta
real(8) function sdelta(stype, x)
Definition
sdelta.f90:10
doccupy.f90
Generated by
1.9.8