The Elk Code
 
Loading...
Searching...
No Matches
gndsteph.f90
Go to the documentation of this file.
1
2! Copyright (C) 2019 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3! E. K. U. Gross. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine gndsteph
7use modmain
8use modphonon
9use modbog
10use modmpi
11use modomp
12implicit none
13! local variables
14integer nmix,nwork
15real(8) dv
16character(64) str
17! allocatable arrays
18real(8), allocatable :: work(:)
19! initialise universal variables
20call init0
21call init1
22call init2
23call readstate
24call genvsig
25call gencore
26call linengy
27call genapwlofr
28call gensocfr
29call genevfsv
30! precise determination of the Fermi energy
32swidth=1.d-5
33call occupy
35! initialise electron-phonon variables
36call initeph
37! size of mixing vector for electron and phonon density matrices (complex array)
38nmix=2*size(duvwx)
39! determine the size of the mixer work array
40nwork=-1
41call mixerifc(mixtype,nmix,duvwx,dv,nwork,work)
42allocate(work(nwork))
43! initialise the mixer
44iscl=0
45call mixerifc(mixtype,nmix,duvwx,dv,nwork,work)
46! set the stop signal to .false.
47tstop=.false.
48! set last self-consistent loop flag
49tlast=.false.
50! only the MPI master process should write files
51if (mp_mpi) then
52! open EPH_INFO.OUT file
53 open(60,file='EPH_INFO.OUT',form='FORMATTED')
54 call writebox(60,"Self-consistent loop started")
55! open EPHGAP.OUT
56 open(64,file='EPHGAP.OUT',form='FORMATTED')
57! open RMSDVS.OUT
58 open(65,file='RMSDVS.OUT',form='FORMATTED')
59! open FACE.OUT
60 open(67,file='FACE.OUT',form='FORMATTED')
61end if
62if (mp_mpi) write(*,*)
63! begin the self-consistent loop
64do iscl=1,maxscl
65 if (mp_mpi) then
66 write(str,'("Loop number : ",I0)') iscl
67 call writebox(60,trim(str))
68 flush(60)
69 write(*,'("Info(gndsteph): self-consistent loop number : ",I4)') iscl
70 end if
71 if (iscl >= maxscl) then
72 if (mp_mpi) then
73 write(60,*)
74 write(60,'("Reached self-consistent loops maximum")')
75 end if
76 tlast=.true.
77 end if
78! determine change in electron and phonon energies
79 call dengyeph
80! solve the electron and phonon eigenvalue equations
81 call eveqneph
82! update the Fermi energy
83 call occupyuv
84 if (mp_mpi) then
85! write the electronic eigenvalues to file
86 call writeevaluv
87! write the phononic eigenvalues to file
88 call writeevalwx
89! write the Fermi energy to file
90 call writefermi
91 write(60,*)
92 write(60,'("Energies :")')
93 write(60,'(" Fermi",T30,": ",G24.14)') efermi
94 write(60,'(" electronic change",T30,": ",G24.14)') dengye
95 write(60,'(" phononic change",T30,": ",G24.14)') dengyph
96 write(60,'(" sum of changes",T30,": ",G24.14)') dengy
97 write(60,*)
98 write(60,'("Estimated indirect band gap : ",G18.10)') bandgap(1)
99 write(60,'(" from k-point ",I6," to k-point ",I6)') ikgap(1),ikgap(2)
100 write(60,'("Estimated direct band gap : ",G18.10)') bandgap(2)
101 write(60,'(" at k-point ",I6)') ikgap(3)
102 write(60,*)
103 write(60,'("Fermionic anomalous correlation entropy : ",G18.10)') face
104 write(60,*)
105 write(60,'("Electron-phonon scaling factor : ",G18.10)') ephscf(1)
106! write estimated indirect band gap
107 write(64,'(G24.14)') bandgap(1)
108 flush(64)
109! write the fermionic anomalous correlation entropy
110 write(67,'(G18.10)') face
111 flush(67)
112 end if
113! mix the old and new electron and phonon density matrices
114 call mixerifc(mixtype,nmix,duvwx,dv,nwork,work)
115! adjust the electron-phonon term scale factor towards 1
116 ephscf(1)=(1.d0-ephscf(2))*ephscf(1)+ephscf(2)
117! exit self-consistent loop if required
118 if (tlast) goto 10
119! check for convergence
120 if (iscl >= 2) then
121 if (mp_mpi) then
122 write(60,*)
123 write(60,'("RMS change in density matrices (target) : ",G18.10," (",&
124 &G18.10,")")') dv,epspot
125 write(65,'(G18.10)') dv
126 flush(65)
127 if (dv < epspot) then
128 write(60,*)
129 write(60,'("Convergence targets achieved")')
130 tlast=.true.
131 end if
132 end if
133 end if
134! check for STOP file
135 call checkstop
136 if (tstop) tlast=.true.
137! broadcast tlast from master process to all other processes
138 call mpi_bcast(tlast,1,mpi_logical,0,mpicom,ierror)
139! reset the OpenMP thread variables
140 call omp_reset
141end do
14210 continue
143if (mp_mpi) then
144 call writebox(60,"Self-consistent loop stopped")
145! close the EPH_INFO.OUT file
146 close(60)
147! close the EPHGAP.OUT file
148 close(64)
149! close the RMSDVS.OUT file
150 close(65)
151! close the FACE.OUT file
152 close(67)
153end if
154deallocate(work)
155! synchronise MPI processes
156call mpi_barrier(mpicom,ierror)
157end subroutine
158
subroutine checkstop
Definition checkstop.f90:7
subroutine dengyeph
Definition dengyeph.f90:7
subroutine eveqneph
Definition eveqneph.f90:7
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine gencore
Definition gencore.f90:10
subroutine genevfsv
Definition genevfsv.f90:7
subroutine gensocfr
Definition gensocfr.f90:7
subroutine genvsig
Definition genvsig.f90:10
subroutine gndsteph
Definition gndsteph.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine initeph
Definition initeph.f90:7
subroutine linengy
Definition linengy.f90:10
subroutine mixerifc(mtype, n, v, dv, nwork, work)
Definition mixerifc.f90:7
complex(8), dimension(:), allocatable, target duvwx
Definition modbog.f90:9
real(8) face
Definition modbog.f90:27
real(8) epspot
Definition modmain.f90:1058
real(8) swidth0
Definition modmain.f90:892
real(8) efermi
Definition modmain.f90:904
integer, dimension(3) ikgap
Definition modmain.f90:914
real(8), dimension(2) bandgap
Definition modmain.f90:912
integer mixtype
Definition modmain.f90:695
real(8) swidth
Definition modmain.f90:892
logical tstop
Definition modmain.f90:1052
logical tlast
Definition modmain.f90:1050
integer maxscl
Definition modmain.f90:1046
integer iscl
Definition modmain.f90:1048
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
subroutine omp_reset
Definition modomp.f90:71
real(8), dimension(2) ephscf
real(8) dengye
real(8) dengy
real(8) dengyph
subroutine occupy
Definition occupy.f90:10
subroutine occupyuv
Definition occupyuv.f90:7
subroutine readstate
Definition readstate.f90:10
subroutine writebox(fnum, str)
Definition writebox.f90:7
subroutine writeevaluv
subroutine writeevalwx
subroutine writefermi