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
6subroutine doccupy
7use modmain
8use modphonon
9implicit none
10! local variables
11integer, parameter :: maxit=1000
12integer ik,jk,ist,it
13real(8) de0,de1,de
14real(8) dchg,x,dx,t1
15! external functions
16real(8), external :: sdelta
17if (.not.tphq0) return
18de0=1.d6
19de1=-1.d6
20do 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
26end do
27t1=1.d0/swidth
28do 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
46end do
47write(*,*)
48write(*,'("Warning(doccupy): could not find Fermi energy derivative")')
49write(*,'(" for s.c. loop ",I5)') iscl
50end subroutine
51
subroutine doccupy
Definition doccupy.f90:7
real(8) efermi
Definition modmain.f90:904
integer nkptnr
Definition modmain.f90:463
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
real(8) swidth
Definition modmain.f90:892
integer nstsv
Definition modmain.f90:886
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer stype
Definition modmain.f90:888
integer iscl
Definition modmain.f90:1048
real(8) wkptnr
Definition modmain.f90:477
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
real(8), dimension(:,:), allocatable doccsv
logical tphq0
Definition modphonon.f90:17
real(8), dimension(:,:), allocatable devalsv
real(8) defermi
real(8) function sdelta(stype, x)
Definition sdelta.f90:10