-
Notifications
You must be signed in to change notification settings - Fork 16
/
Copy pathtotal_energy.f90
287 lines (238 loc) · 9.57 KB
/
total_energy.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
#include "alias.inc"
module total_energy
use parameters, only: energy, incar, kpoints, poscar
use mpi_setup
use sorting, only : get_sort_index_1D
use do_math, only : idx2Di, idx2Dj
use print_io
contains
subroutine get_total_energy(ETBA, PINPT, PKPTS, PGEOM, flag_opt_ef, flag_report, E_F_)
implicit none
type(energy) :: ETBA
type(incar) :: PINPT
type(poscar) :: PGEOM
type(kpoints) :: PKPTS
integer*4 nkp, nband, nspin
integer*4 mpierr
real*8 eltemp
real*8 nelect_ref
real*8 degen ! degeneracy
real*8 E_F, E0
real*8 OCC(PGEOM%nband*PINPT%nspin,PKPTS%nkpoint)
logical flag_quit, flag_opt_ef, flag_report
real*8,intent(inout),optional :: E_F_
if(.not. PINPT%flag_get_total_energy) return
if(flag_report) then
write(message,'(A)')' ' ; write_msg
write(message,'(A)')' #--- START: TOTAL ENERGY EVALUATION -----------' ; write_msg
endif
! initialize
ETBA%E_BAND = 0d0
ETBA%E_TOT = 0d0
degen = 1d0
ETBA%E_F = 0d0 ! temporary
! if present E_F_, use it and return it (if flag_opt_ef -> update) with E_F_,
! otherwise, use ETBA%E_F and (if flag_opt_ef -> update) save it to ETBA%E_F
if(present(E_F_)) then
E_F = E_F_
else
E_F = ETBA%E_F
endif
eltemp = PINPT%electronic_temperature ! default = 0k
nkp = PKPTS%nkpoint
nband = PGEOM%nband
nspin = PINPT%nspin
nelect_ref = sum(PGEOM%nelect) ! reference number of electrons
if(PINPT%ispin .eq. 1) degen = 2d0
! optimize Fermi level with given number of electrons nelect_ref
if(flag_opt_ef) then
call get_fermi_level(ETBA%E, E_F, nelect_ref, nband, nspin, nkp, eltemp, degen)
if(present(E_F_)) then
E_F_ = E_F
else
ETBA%E_F = E_F
endif
endif
! evaluate occupation with Fermi level E_F, F_OCC
call get_occupation(ETBA%E, E_F, ETBA%F_OCC, nband, nspin, nkp, degen, eltemp)
! evaluate band energy based on F_OCC
ETBA%E_BAND = sum(band_energy(ETBA%E, ETBA%F_OCC))/nkp
! evaluate band energy based on F_OCC with eltemp = 0d0
call get_occupation(ETBA%E, E_F, OCC, nband, nspin, nkp, degen, 0d0)
E0 = sum(band_energy(ETBA%E, OCC))/nkp
if(flag_report) then
write(message,'(A,F20.8, A)') ' -Fermi level E_F (in eV) : ', E_F ; write_msg
write(message,'(A,F20.8, A)') ' -Number of electrons (with E_F) : ', sum(ETBA%F_OCC)/nkp ; write_msg
write(message,'(A,F20.8, A)') ' -Number of electrons (NELECT) : ', nelect_ref ; write_msg
write(message,'(A)') ' -Total energy components (in eV) ' ; write_msg
write(message,'(A,F10.4, A,F20.8)')' *- Band energy (ELTEMP = ',eltemp,') : ', ETBA%E_BAND ; write_msg
write(message,'(A,F10.4, A,F20.8)')' *- Band energy (ELTEMP -> ',0.0d0 ,') : ', E0 ; write_msg
write(message,'(A )') ' __________________________________________' ; write_msg
write(message,'(A,F20.8, A)') ' #- Total energy (in eV) : ', ETBA%E_BAND ; write_msg
write(message,'(A)')' #--- END: TOTAL ENERGY EVALUATION -----------' ; write_msg
write(message,'(A)')' ' ; write_msg
endif
return
endsubroutine
subroutine get_fermi_level(E, E_F, nelect_ref, nband, nspin, nkp, eltemp, degen)
implicit none
integer*4 ii, ik, ie
integer*4 nkp, nband, nspin
integer*4 mpierr
real*8 eltemp, de_cut
real*8 nelect
real*8 nelect_ref
real*8 grad_nelect
real*8 degen ! degeneracy
real*8 E(nband*nspin,nkp)
real*8 E_F
logical flag_quit
nelect = 0d0
flag_quit = .false.
grad_nelect = 0d0
de_cut = 0.0000001 ! cutoff for the convergence of Fermi level evaluation
!step0. Rough estimation of Fermi level with mid-gap position
call get_midgap(E, E_F, nelect_ref, nband, nspin, nkp, degen)
!step1. get nelect with initial E_F = 0d0
call get_nelect(E, E_F, nelect, nband, nspin, nkp, eltemp, degen)
if( abs(nelect - nelect_ref) .lt. de_cut) then
flag_quit = .true.
else
flag_quit = .false.
call get_grad_nelect(E, E_F, grad_nelect, nband, nspin, nkp, eltemp, degen)
E_F = E_F - (nelect - nelect_ref)/grad_nelect
call get_nelect(E, E_F, nelect, nband, nspin, nkp, eltemp, degen)
endif
!step2. if nelect = nelec_ref -> quit, otherwise, recalculate with updated E_F
! based on gradient of nelect w.r.t chemical potential
ii = 0
do while (.not. flag_quit)
ii = ii + 1
if( abs(nelect - nelect_ref) .lt. de_cut ) then
flag_quit = .true.
else
flag_quit = .false.
call get_grad_nelect(E, E_F, grad_nelect, nband, nspin, nkp, eltemp, degen)
E_F = E_F - (nelect - nelect_ref)/grad_nelect
call get_nelect(E, E_F, nelect, nband, nspin, nkp, eltemp, degen)
endif
enddo
return
endsubroutine
elemental real*8 function band_energy(E, occ)
implicit none
real*8, intent(in) :: E, occ
band_energy = E * occ
return
endfunction
subroutine get_occupation(E, E_F, occ, nband, nspin, nkp, degen, eltemp)
use do_math, only : fdirac
implicit none
integer*4 nband, nspin, nkp
integer*4 ik, ie, is
integer*4 mpierr
real*8 E_F
real*8 eltemp, degen
real*8 E(nband*nspin, nkp)
real*8 occ(nband * nspin, nkp), occ_(nband * nspin, nkp)
integer*4 ourjob(nprocs), ourjob_disp(0:nprocs-1)
occ_ = 0d0
call mpi_job_distribution_chain(nkp, nprocs, ourjob, ourjob_disp)
do ik = sum(ourjob(1:myid))+1, sum(ourjob(1:myid+1))
do is = 1, nspin
occ_(1+nband*(is-1):nband*is,ik) = fdirac( E(1+nband*(is-1):nband*is,ik), E_F, eltemp ) * degen
enddo
enddo
#ifdef MPI
call MPI_ALLREDUCE(occ_, occ, size(occ), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
#else
occ = occ_
#endif
return
endsubroutine
subroutine get_nelect(E, E_F, nelect, nband, nspin, nkp, eltemp, degen)
use do_math, only : fdirac
implicit none
integer*4 nband, nspin, nkp
integer*4 ik, ie, is
integer*4 mpierr
real*8 dos_tot
real*8 eltemp, degen
real*8 E(nband*nspin, nkp)
real*8 nelect
real*8 E_F
integer*4 ourjob(nprocs), ourjob_disp(0:nprocs-1)
call mpi_job_distribution_chain(nkp, nprocs, ourjob, ourjob_disp)
dos_tot = 0d0
do ik = sum(ourjob(1:myid))+1, sum(ourjob(1:myid+1))
do is = 1, nspin
dos_tot = dos_tot + sum( fdirac( E(1+nband*(is-1):nband*is,ik), E_F, eltemp ) ) * degen
enddo
enddo
#ifdef MPI
call MPI_ALLREDUCE( dos_tot/real(nkp), nelect, 1, MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
#else
nelect = dos_tot/real(nkp)
#endif
return
endsubroutine
subroutine get_grad_nelect(E, E_F, grad_nelect, nband, nspin, nkp, eltemp, degen)
use do_math, only : grad_u_fdirac
implicit none
integer*4 nband, nspin, nkp
integer*4 ik, ie, is
integer*4 mpierr
real*8 eltemp, degen
real*8 E(nband*nspin, nkp)
real*8 grad_nelect, grad_nelect_
real*8 E_F
integer*4 ourjob(nprocs), ourjob_disp(0:nprocs-1)
call mpi_job_distribution_chain(nkp, nprocs, ourjob, ourjob_disp)
grad_nelect = 0d0
do ik = sum(ourjob(1:myid))+1, sum(ourjob(1:myid+1))
do is = 1, nspin
grad_nelect = grad_nelect + sum( grad_u_fdirac( E(1+nband*(is-1):nband*is,ik), E_F, eltemp ) ) * degen
enddo
enddo
#ifdef MPI
call MPI_ALLREDUCE(grad_nelect, grad_nelect_, 1, MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr)
grad_nelect = grad_nelect_ / real(nkp)
#else
grad_nelect = grad_nelect / real(nkp)
#endif
return
endsubroutine
subroutine get_midgap(E, E_F, nelect_ref, nband, nspin, nkp, degen)
use sorting, only : get_sort_index_1D
use do_math, only : idx2Di, idx2Dj
implicit none
integer*4 ii, ie, ik
integer*4 nband, nspin, nkp, neig_tot
integer*4 idx_eig1D(nband*nspin*nkp)
real*8, intent(in) :: E(nband*nspin, nkp)
real*8 nelect_ref
real*8 nelect
real*8 E1, E2
real*8 E_F
real*8 degen
neig_tot = nband*nspin*nkp
nelect = 0d0
ii = 0
! sort eigenvalue in ascending order
call get_sort_index_1D(idx_eig1D, reshape(E,[neig_tot]), neig_tot, 'ascending ')
do while (nelect .lt. nelect_ref)
ii = ii + 1
nelect = nelect + 1d0/real(nkp) * degen
ie = idx2Di(idx_eig1D(ii),nband*nspin)
ik = idx2Dj(idx_eig1D(ii),nband*nspin)
enddo
! VBM
E1 = E(ie, ik)
! CBM
ie = idx2Di(idx_eig1D(ii),nband*nspin)
ik = idx2Dj(idx_eig1D(ii),nband*nspin)
E2 = E(ie, ik)
E_F = (E1 + E2)*0.5d0
return
endsubroutine
endmodule