The Elk Code
 
Loading...
Searching...
No Matches
gwefermi.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 P. Elliott, 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 gwefermi
7use modmain
8use modgw
9use modmpi
10use modomp
11implicit none
12! local variables
13logical done
14integer, parameter :: maxit=1000
15integer ik,ist,it,nthd
16real(8) e0,e1,e,chg,chgk
17if (mp_mpi) then
18 write(*,*)
19 write(*,'("Info(gwefermi): finding the GW Fermi energy")')
20end if
21! find minimum and maximum eigenvalues
22e0=evalsv(1,1)
23e1=e0
24do ik=1,nkpt
25 do ist=1,nstsv
26 e=evalsv(ist,ik)
27 if (e < e0) e0=e
28 if (e > e1) e1=e
29 end do
30end do
31done=.false.
32do it=1,maxit
33 if (mp_mpi.and.(mod(it,10) == 0)) then
34 write(*,'("Info(gwefermi): done ",I4," iterations")') it
35 end if
36 efermi=0.5d0*(e0+e1)
37 chg=0.d0
38! begin parallel loop over k-points
39 call holdthd(nkpt/np_mpi,nthd)
40!$OMP PARALLEL DO DEFAULT(SHARED) &
41!$OMP PRIVATE(chgk) REDUCTION(+:chg) &
42!$OMP SCHEDULE(DYNAMIC) &
43!$OMP NUM_THREADS(nthd)
44 do ik=1,nkpt
45! distribute among MPI processes
46 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
47 call gwchgk(ik,chgk)
48 chg=chg+chgk
49 end do
50!$OMP END PARALLEL DO
51 call freethd(nthd)
52! add charge from each process and redistribute
53 if (np_mpi > 1) then
54 call mpi_allreduce(mpi_in_place,chg,1,mpi_double_precision,mpi_sum,mpicom, &
55 ierror)
56 end if
57 if (chg < chgval) then
58 e0=efermi
59 else
60 e1=efermi
61 end if
62! check for convergence
63 if ((e1-e0) < 1.d-12) done=.true.
64! broadcast done from master process to all other processes
65 call mpi_bcast(done,1,mpi_logical,0,mpicom,ierror)
66 if (done) return
67end do
68write(*,*)
69write(*,'("Warning(gwefermi): could not find GW Fermi energy")')
70end subroutine
71
subroutine gwchgk(ik, chgk)
Definition gwchgk.f90:7
subroutine gwefermi
Definition gwefermi.f90:7
Definition modgw.f90:6
real(8) efermi
Definition modmain.f90:904
real(8) chgval
Definition modmain.f90:722
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
logical mp_mpi
Definition modmpi.f90:17
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106