Octopus
smear.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19! Some pieces glued from PW/w(0,1)gauss.f90 from PWSCF
20
21#include "global.h"
22
23module smear_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
29 use math_oct_m, only: factorial
32 use parser_oct_m
35 use sort_oct_m
36 use unit_oct_m
39
40 implicit none
41
42 private
43 public :: &
44 smear_t, &
45 smear_init, &
46 smear_copy, &
55
56 type smear_t
57 private
58 integer, public :: method
59 real(real64), public :: dsmear
60 real(real64) :: dsmear_cond
61 real(real64), public :: e_fermi
62 real(real64) :: e_fermi_cond
63
64 integer, public :: el_per_state
65 real(real64), public :: ef_occ
66 real(real64) :: ef_occ_cond
67 logical, public :: integral_occs
68 integer, public :: MP_n
69 integer :: fermi_count
70 integer :: nik_factor
71 integer :: nspins
72
73 type(simplex_t), pointer :: tetra_simplex => null()
74
75 ! PhotonDoping
76 logical, public:: photodop
77 real(real64) :: nephotodop
78 integer :: photodop_bandmin
79
80 end type smear_t
81
82 integer, parameter, public :: &
83 SMEAR_SEMICONDUCTOR = 1, &
85 smear_cold = 3, &
87 smear_spline = 5, &
88 smear_fixed_occ = 6, &
89 smear_gaussian = 7, &
90 smear_lorentzian = 8, &
91 smear_tetrahedra = 21, &
93
94 real(real64), parameter :: TOL_SMEAR = 1e-6_real64
95
96contains
97
98 !--------------------------------------------------
99 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
100 type(smear_t), intent(out) :: this
101 type(namespace_t), intent(in) :: namespace
102 integer, intent(in) :: ispin
103 logical, intent(in) :: fixed_occ
104 logical, intent(in) :: integral_occs
105 type(kpoints_t), intent(in) :: kpoints
106
107 push_sub(smear_init)
108
109 this%integral_occs = integral_occs
110
111 !! The numbers 21 and 22 for tetrahedra options are reused in src/electrons/dos.F90
112 !%Variable SmearingFunction
113 !%Type integer
114 !%Default semiconducting
115 !%Section States
116 !%Description
117 !% This is the function used to smear the electronic occupations.
118 !% It is ignored if the <tt>Occupations</tt> block is set.
119 !%Option semiconducting 1
120 !% Semiconducting occupations, <i>i.e.</i> the lowest lying states are occupied
121 !% until no more electrons are left.
122 !%Option fermi_dirac 2
123 !% Simple Fermi-Dirac distribution. In this case, <tt>Smearing</tt> has
124 !% the meaning of an electronic temperature. DN Mermin, <i>Phys. Rev.</i> <b>137</b>, A1441 (1965).
125 !%Option cold_smearing 3
126 !% N Marzari, D Vanderbilt, A De Vita, and MC Payne, <i>Phys. Rev. Lett.</i> <b>82</b>, 3296 (1999).
127 !%Option methfessel_paxton 4
128 !% M Methfessel and AT Paxton, <i>Phys. Rev. B</i> <b>40</b>, 3616 (1989).
129 !% In this case, the variable <tt>SmearingMPOrder</tt> sets the order of the smearing.
130 !% Occupations may be negative.
131 !%Option spline_smearing 5
132 !% Nearly identical to Gaussian smearing.
133 !% JM Holender, MJ Gillan, MC Payne, and AD Simpson, <i>Phys. Rev. B</i> <b>52</b>, 967 (1995).
134 !%Option gaussian_smearing 7
135 !% Gaussian smearing.
136 !%Option lorentzian_smearing 8
137 !% Lorentzian smearing.
138 !%Option tetrahedra 21
139 !% Linear tetrahedron method as in P. E. Bloechl, et al., <i>Phys. Rev. B</i> <b>49</b>, 16223 (1994).
140 !% Requires a regular Monkhorst-Pack generated by Octopus.
141 !%Option tetrahedra_opt 22
142 !% Improved tetrahedron method as in M. Kawamura, et al., <i>Phys. Rev. B</i> <b>89</b>, 094515 (2014).
143 !% Requires a regular Monkhorst-Pack generated by Octopus.
144 !%End
145 if (fixed_occ) then
146 this%method = smear_fixed_occ
147 else
148 call parse_variable(namespace, 'SmearingFunction', smear_semiconductor, this%method)
149 if (.not. varinfo_valid_option('SmearingFunction', this%method)) call messages_input_error(namespace, 'SmearingFunction')
150 call messages_print_var_option('SmearingFunction', this%method, namespace=namespace)
151 end if
152
153 !%Variable Smearing
154 !%Type float
155 !%Default 0.1 eV
156 !%Section States
157 !%Description
158 !% If <tt>Occupations</tt> is not set, <tt>Smearing</tt> is the
159 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons
160 !% among the existing states.
161 !%End
162 this%dsmear = 1e-14_real64
163 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
164 call parse_variable(namespace, 'Smearing', 0.1_real64 / (m_two * p_ry), this%dsmear, units_inp%energy)
165 end if
167 !%Variable PhotoDopingSmearing
168 !%Type float
169 !%Default 0.1 eV
170 !%Section States
171 !%Description
172 !% If <tt>Occupations</tt> is not set, <tt>PhotoDopingSmearing</tt> is the
173 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons in conduction bands
174 !% among the existing states.
175 !%End
176 this%dsmear_cond = this%dsmear
177 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
178 call parse_variable(namespace, 'PhotoDopingSmearing', 0.1_real64 / (m_two * p_ry), this%dsmear_cond, units_inp%energy)
179 end if
180
181 !%Variable PhotoDopingNumElectrons
182 !%Type float
183 !%Default 0
184 !%Section States
185 !%Description
186 !% This variable allows to set the number of valence electrons taken out to put
187 !% into the conduction bands. The method follows Marini and Calandra, PRB 104, 144103 (2021).
188 !%End
189 call parse_variable(namespace, 'PhotoDopingNumElectrons', m_zero, this%nephotodop)
190
191 if (this%nephotodop > m_min_occ) then
192 this%photodop=.true.
193 else
194 this%photodop=.false.
195 endif
196
197
198 !%Variable PhotoDopingBand
199 !%Type integer
200 !%Default 1
201 !%Section States
202 !%Description
203 !% This variable specifies where the minium (starting) conduction band index for the photodoping electrons.
204 !%End
205 call parse_variable(namespace, 'PhotoDopingBand', 1, this%photodop_bandmin)
206
207
208 call messages_obsolete_variable(namespace, 'ElectronicTemperature', 'Smearing')
209
210 if (ispin == 1) then ! unpolarized
211 this%el_per_state = 2
212 else
213 this%el_per_state = 1
214 end if
215
216 if (ispin == 2) then
217 this%nspins = 2
218 else
219 this%nspins = 1
220 end if
221
222 if (this%method == smear_semiconductor) then
223 this%nik_factor = kpoints_kweight_denominator(kpoints)
224
225 if (this%nik_factor == 0) then
226 message(1) = "k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
227 call messages_fatal(1, namespace=namespace)
228 end if
229 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
230 if (.not. associated(kpoints%simplex)) then
231 message(1) = "Tetrahedron smearing requires a simplex mesh for the k-point grid."
232 call messages_fatal(1, namespace=namespace)
233 end if
234
235 this%tetra_simplex => kpoints%simplex
236 end if
237
238 this%MP_n = 0
239 if (this%method == smear_methfessel_paxton) then
240 !%Variable SmearingMPOrder
241 !%Type integer
242 !%Default 1
243 !%Section States
244 !%Description
245 !% Sets the order of the Methfessel-Paxton smearing function.
246 !%End
247 call parse_variable(namespace, 'SmearingMPOrder', 1, this%MP_n)
248 end if
249
250 pop_sub(smear_init)
251 end subroutine smear_init
252
253
254 !--------------------------------------------------
255 subroutine smear_copy(to, from)
256 type(smear_t), intent(out) :: to
257 type(smear_t), intent(in) :: from
258
259 push_sub(smear_copy)
260
261 to%method = from%method
262 to%dsmear = from%dsmear
263 to%e_fermi = from%e_fermi
264 to%el_per_state = from%el_per_state
265 to%ef_occ = from%ef_occ
266 to%MP_n = from%MP_n
267 to%fermi_count = from%fermi_count
268 to%nik_factor = from%nik_factor
269 to%tetra_simplex => from%tetra_simplex
270
271 pop_sub(smear_copy)
272 end subroutine smear_copy
273
274
275 !--------------------------------------------------
276 subroutine smear_find_fermi_energy(this, namespace, eigenvalues, occupations, &
277 qtot, nik, nst, kweights)
278 type(smear_t), intent(inout) :: this
279 type(namespace_t), intent(in) :: namespace
280 real(real64), intent(in) :: eigenvalues(:,:), occupations(:,:)
281 real(real64), intent(in) :: qtot, kweights(:)
282 integer, intent(in) :: nik, nst
283
284 real(real64), parameter :: tol = 1.0e-10_real64
285 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
286 real(real64) :: sumq_frac
287 logical :: conv
288 real(real64), allocatable :: eigenval_list(:)
289 integer, allocatable :: k_list(:), reorder(:)
290 integer :: fermi_count_up, fermi_count_down
291
293
294 maxq = this%el_per_state * nst * this%nspins
295 if (maxq - qtot <= -tol) then ! not enough states
296 message(1) = 'Not enough states'
297 write(message(2),'(6x,a,f12.6,a,i10)')'(total charge = ', qtot, &
298 ' max charge = ', maxq
299 call messages_fatal(2, namespace=namespace)
300 end if
301
302 conv = .true.
303 if (this%method == smear_fixed_occ) then ! Fermi energy: last of occupied states
304 ist_cycle: do ist = nst, 1, -1
305 if (any(occupations(ist, :) > m_min_occ)) then
306 this%e_fermi = eigenvalues(ist, 1)
307 this%ef_occ = occupations(ist, 1) / this%el_per_state
308 do ik = 2, nik
309 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) > m_min_occ) then
310 this%e_fermi = eigenvalues(ist, ik)
311 this%ef_occ = occupations(ist, ik) / this%el_per_state
312 end if
313 end do
314 exit ist_cycle
315 end if
316 end do ist_cycle
317
318 else if (this%method == smear_semiconductor) then
319 ! first we sort the eigenvalues
320 safe_allocate(eigenval_list(1:nst * nik))
321 safe_allocate( k_list(1:nst * nik))
322 safe_allocate( reorder(1:nst * nik))
323
324 iter = 1
325 do ist = 1, nst
326 do ik = 1, nik
327 eigenval_list(iter) = eigenvalues(ist, ik)
328 k_list(iter) = ik
329 reorder(iter) = iter
330 iter = iter + 1
331 end do
332 end do
333
334 call sort(eigenval_list, reorder)
335
336 sumq_int = floor(qtot) * this%nik_factor
337 sumq_frac = qtot * this%nik_factor - sumq_int
338 if ( sumq_frac + tol > this%nik_factor ) then
339 sumq_int = sumq_int + this%nik_factor
340 sumq_frac = sumq_frac - this%nik_factor
341 end if
342 if (sumq_frac < tol) sumq_frac = m_zero
343
344 do iter = 1, nst * nik
345 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor + m_half)
346 if (.not. weight > 0) cycle
347 this%e_fermi = eigenval_list(iter)
348 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
349
350 if (sumq_int - weight * this%el_per_state <= 0) then
351 ! count how many occupied states are at the fermi level,
352 ! this is required later to fill the states
353 this%fermi_count = 1
354 fermi_count_up = 1
355 fermi_count_down = 1
356 sum_weight = weight
357 do
358 if (iter - fermi_count_down < 1) exit
359 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear) exit
360 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor + m_half)
361 if (weight > m_epsilon) then
362 sumq_int = sumq_int + weight * this%el_per_state
363 sum_weight = sum_weight + weight
364 end if
365 fermi_count_down = fermi_count_down + 1
366 end do
367 do
368 if (iter + fermi_count_up > nst*nik) exit
369 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear) exit
370 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor + m_half)
371 if (weight > m_epsilon) then
372 sum_weight = sum_weight + weight
373 end if
374 fermi_count_up = fermi_count_up + 1
375 end do
376 this%fermi_count = fermi_count_up + fermi_count_down - 1
377 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
378 exit
379 end if
380
381 sumq_int = sumq_int - weight * this%el_per_state
382 end do
383 assert(this%ef_occ < m_one + 10.0_real64*m_epsilon)
384
385 safe_deallocate_a(eigenval_list)
386 safe_deallocate_a(k_list)
387 safe_deallocate_a(reorder)
388
389 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
390 if (this%photodop) then
391 call messages_not_implemented('Photodoping for SmearingFunction = tetrahedra/tetrahedra_opt')
392 end if
393 call bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, qtot, this%e_fermi)
394 else ! bisection
395 if (this%photodop) then
396 ! find fermi energy for valence electrons
397 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
398 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
399 this%e_fermi, this%ef_occ)
400
401 ! find fermi energy for conduction electrons
402 call bisection_find_fermi_energy(this, namespace, this%dsmear_cond, tol, &
403 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
404 this%e_fermi_cond, this%ef_occ_cond)
405 else
406 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
407 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
408 end if
409 end if
410
412 end subroutine smear_find_fermi_energy
413
414 ! ---------------------------------------------------------
415 subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, &
416 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
417 type(smear_t), intent(inout) :: this
418 type(namespace_t), intent(in) :: namespace
419 real(real64), intent(in) :: dsmear_in
420 real(real64), intent(in) :: tol
421 integer, intent(in) :: nik
422 real(real64), intent(in) :: eigenvalues(:,:)
423 real(real64), intent(in) :: kweights(:), q_in
424 integer, intent(in) :: start_band, end_band
425 real(real64), intent(out) :: e_fermi, ef_occ
426
427 integer, parameter :: nitmax = 200
428 integer :: ist, ik, iter
429 real(real64) :: drange, dsmear, emin, emax, xx, sumq
430 logical :: conv
431
433
434 dsmear = max(1e-14_real64, dsmear_in)
435 drange = dsmear * sqrt(-log(tol * 0.01_real64))
436
437 emin = minval(eigenvalues) - drange
438 emax = maxval(eigenvalues) + drange
439
440 do iter = 1, nitmax
441 e_fermi = m_half * (emin + emax)
442 sumq = m_zero
443
444 do ik = 1, nik
445 do ist = start_band, end_band
446 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
447 sumq = sumq + kweights(ik) * this%el_per_state * &
448 smear_step_function(this, xx)
449 end do
450 end do
451
452 conv = (abs(sumq - q_in) <= tol)
453 if (conv) exit
454
455 if (sumq <= q_in) emin = e_fermi
456 if (sumq >= q_in) emax = e_fermi
457
458 ef_occ = smear_step_function(this, m_zero)
459 end do
460
461 if (.not. conv) then
462 message(1) = 'Fermi: did not converge.'
463 call messages_fatal(1, namespace=namespace)
464 end if
465
467 end subroutine bisection_find_fermi_energy
468
469 ! ---------------------------------------------------------
489 subroutine bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, q_in, e_fermi)
490 type(smear_t), intent(inout) :: this
491 type(namespace_t), intent(in) :: namespace
492 real(real64), intent(in) :: tol
493 integer, intent(in) :: nik, nst
494 real(real64), intent(in) :: eigenvalues(:,:)
495 real(real64), intent(in) :: q_in
496 real(real64), intent(out) :: e_fermi
497
498 integer, parameter :: nitmax = 200
499 integer :: iter
500 real(real64) :: emin, emax, sumq
501 logical :: conv
502
503 integer :: ist, ii, ll, is, ik, ns, nik_red, nvertices
504 real(real64) :: e_simplex(20), w_simplex(20), dos_simplex, simplex_unit_volume
505
506 push_sub_with_profile(bisection_find_fermi_energy_tetra)
507
508 assert(associated(this%tetra_simplex))
509 assert(this%tetra_simplex%n_simplices > 0)
510 assert(this%tetra_simplex%sdim > 0)
511
512 ns = this%nspins
513 assert(ns > 0)
514 assert(mod(nik, ns) == 0)
515
516 nik_red = nik / ns
517 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
518
519 nvertices = this%tetra_simplex%rdim + 1
520 assert(nvertices > 0)
521 assert(nvertices <= this%tetra_simplex%sdim)
522
523 simplex_unit_volume = factorial(this%tetra_simplex%rdim)
524
525 emin = minval(eigenvalues)
526 emax = maxval(eigenvalues)
527 if (.not. (emin < emax)) then
528 write(message(1), '(a)') 'Tetra bisection: invalid energy bracket from eigenvalues.'
529 call messages_fatal(1, namespace=namespace)
530 end if
531
532 do iter = 1, nitmax
533 e_fermi = m_half * (emin + emax)
534 sumq = m_zero
535 !$omp parallel do collapse(2) private(ist, ii, is, ll, ik, e_simplex, w_simplex, dos_simplex) reduction(+:sumq)
536 do ist = 1, nst
537 do ii = 1, this%tetra_simplex%n_simplices
538 do is = 0, ns - 1
539 do ll = 1, this%tetra_simplex%sdim
540 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
541 e_simplex(ll) = eigenvalues(ist, ik + is)
542 end do
543 w_simplex(:) = m_zero
544
545 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), e_fermi, &
546 w_simplex(1:this%tetra_simplex%sdim), dos_simplex)
547
548 do ll = 1, nvertices
549 sumq = sumq + this%el_per_state * simplex_unit_volume * w_simplex(ll) / this%tetra_simplex%n_simplices
550 end do
551 end do
552 end do
553 end do
554 !$omp end parallel do
555
556 conv = (abs(sumq - q_in) <= tol)
557 if (conv) exit
558
559 if (sumq <= q_in) emin = e_fermi
560 if (sumq >= q_in) emax = e_fermi
561 end do
562
563 if (.not. conv) then
564 message(1) = 'Fermi (tetrahedra): did not converge.'
565 call messages_fatal(1, namespace=namespace)
566 end if
567
568 pop_sub_with_profile(bisection_find_fermi_energy_tetra)
570
571 ! ---------------------------------------------------------
572 subroutine smear_fill_occupations(this, eigenvalues, occupations, kweights, nik, nst)
573 type(smear_t), intent(in) :: this
574 real(real64), intent(in) :: eigenvalues(:,:)
575 real(real64), intent(inout) :: occupations(:,:)
576 real(real64), intent(in) :: kweights(:)
577 integer, intent(in) :: nik, nst
578
579 integer :: ik, ist, ifermi
580 real(real64) :: dsmear, xx, dsmear_cond
581
582 push_sub(smear_fill_occupations)
583
584 if (this%method == smear_fixed_occ) then
585 ! do nothing
586 else if (this%method == smear_semiconductor) then
587 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
588
589 ifermi = 0
590 do ik = 1, nik
591 do ist = 1, nst
592 xx = eigenvalues(ist, ik) - this%e_fermi
593 if (xx < -tol_smear) then
594 occupations(ist, ik) = this%el_per_state
595 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count) then
596 occupations(ist, ik) = this%ef_occ * this%el_per_state
597 ifermi = ifermi + 1
598 else
599 occupations(ist, ik) = m_zero
600 end if
601
602 end do
603 end do
604 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
605 if (this%photodop) then
606 call messages_not_implemented('Photodoping for SmearingFunction = tetrahedra/tetrahedra_opt')
607 end if
608 block
609 integer :: ii, ll, is, ns, nik_red, nvertices
610 real(real64) :: e_simplex(20), w_simplex(20), dos_simplex, simplex_unit_volume
611
612 assert(associated(this%tetra_simplex))
613 assert(this%tetra_simplex%n_simplices > 0)
614 assert(this%tetra_simplex%sdim > 0)
615
616 ns = this%nspins
617 assert(ns > 0)
618 assert(mod(nik, ns) == 0)
619
620 nik_red = nik / ns
621 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
622
623 nvertices = this%tetra_simplex%rdim + 1
624 assert(nvertices > 0)
625 assert(nvertices <= this%tetra_simplex%sdim)
626
627 simplex_unit_volume = factorial(this%tetra_simplex%rdim)
628
629 occupations(:, :) = m_zero
630 !$omp parallel do collapse(2) private(ist, ii, is, ll, ik, e_simplex, w_simplex, dos_simplex) reduction(+:occupations)
631 do ist = 1, nst
632 do ii = 1, this%tetra_simplex%n_simplices
633 do is = 0, ns - 1
634 do ll = 1, this%tetra_simplex%sdim
635 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
636 e_simplex(ll) = eigenvalues(ist, ik + is)
637 end do
638 w_simplex(:) = m_zero
639
640 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), this%e_fermi, &
641 w_simplex(1:this%tetra_simplex%sdim), dos_simplex)
642
643 do ll = 1, nvertices
644 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
645 assert(kweights(ik + is) > m_epsilon)
646 occupations(ist, ik + is) = occupations(ist, ik + is) + this%el_per_state * &
647 simplex_unit_volume * w_simplex(ll) / (this%tetra_simplex%n_simplices * kweights(ik + is))
648 end do
649 end do
650 end do
651 end do
652 !$omp end parallel do
653 end block
654 else
655 dsmear = max(1e-14_real64, this%dsmear)
656 if (this%photodop) then
657 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
658 do ik = 1, nik
659 do ist = 1, nst
660 if (ist < this%photodop_bandmin) then
661 ! valence electrons
662 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
663 else
664 ! conduction electrons
665 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
666 end if
667 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
668 end do
669 end do
670 else
671 do ik = 1, nik
672 do ist = 1, nst
673 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
674 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
675 end do
676 end do
677 end if
678 end if
679
681 end subroutine smear_fill_occupations
682
683 !--------------------------------------------------
684 real(real64) function smear_calc_entropy(this, eigenvalues, &
685 nik, nst, kweights, occ) result(entropy)
686 type(smear_t), intent(inout) :: this
687 real(real64), intent(in) :: eigenvalues(:,:)
688 real(real64), intent(in) :: kweights(:)
689 integer, intent(in) :: nik, nst
690 real(real64), intent(in) :: occ(:, :)
691
692 integer :: ist, ik
693 real(real64) :: dsmear, xx, term, ff
694
695 push_sub(smear_calc_entropy)
696
697 entropy = m_zero
698
699 if (this%method == smear_fixed_occ .or. this%method == smear_semiconductor) then
700 ! Fermi-Dirac entropy, not quite the same as will be obtained with true smearing
701 ! RM Wentzcovitch, JL Martins, and PB Allen, Phys. Rev. B 45, 11372 (1992) eqn (5)
702 ! also N Marzari PhD thesis p 117, http://quasiamore.mit.edu/phd/Marzari_PhD.pdf
703 do ik = 1, nik
704 do ist = 1, nst
705 ff = occ(ist, ik) / this%el_per_state
706 if (ff > m_zero .and. ff < m_one) then
707 term = ff * log(ff) + (1 - ff) * log(1 - ff)
708 else ! we have semiconducting smearing, or perverse occupations as in Methfessel-Paxton
709 term = m_zero
710 end if
711 entropy = entropy - kweights(ik) * this%el_per_state * term
712 end do
713 end do
714 else
715 dsmear = max(1e-14_real64, this%dsmear)
716
717 do ik = 1, nik
718 do ist = 1, nst
719 if (eigenvalues(ist, ik) < huge(m_one)) then
720 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
721 entropy = entropy - kweights(ik) * this%el_per_state * &
722 smear_entropy_function(this, xx)
723 end if
724 end do
725 end do
726 end if
727
728 pop_sub(smear_calc_entropy)
729 end function smear_calc_entropy
730
731
732 ! ---------------------------------------------------------
733 real(real64) function smear_delta_function(this, xx) result(deltaf)
734 type(smear_t), intent(in) :: this
735 real(real64), intent(in) :: xx
736
737 real(real64), parameter :: maxarg = 200.0_real64
738 real(real64) :: xp, arg, hd, hp, aa
739 integer :: ii, ni
740
741 ! no PUSH_SUB, called too often
742
743 ! smear_delta_function is not defined for SMEAR_FIXED_OCC
744 assert(this%method /= smear_fixed_occ)
745
746 deltaf = m_zero
747 select case (this%method)
749 if (abs(xx) <= m_epsilon) then
750 deltaf = this%ef_occ
751 end if
752
753 case (smear_fermi_dirac)
754 if (abs(xx) <= 36.0_real64) then
755 deltaf = m_one / (m_two + exp(-xx) + exp(xx))
756 end if
757
758 case (smear_cold)
759 xp = xx - m_one / sqrt(m_two)
760 arg = min(maxarg, xp**2)
761
762 deltaf = exp(-arg) / sqrt(m_pi) * (m_two - sqrt(m_two) * xx)
763
765 arg = min(maxarg, xx**2)
766 deltaf = exp(-arg) / sqrt(m_pi)
767
768 if (this%MP_n > 0) then ! recursion
769 hd = m_zero
770 hp = exp(-arg)
771 ni = 0
772 aa = m_one / sqrt(m_pi)
773 do ii = 1, this%MP_n
774 hd = m_two * xx * hp - m_two * ni * hd
775 ni = ni + 1
776 aa = -aa / (m_four * ii)
777 hp = m_two * xx * hd - m_two * ni * hp
778 ni = ni + 1
779 deltaf = deltaf + aa * hp
780 end do
781 end if
782
783 case (smear_spline)
784 xp = abs(xx) + m_one / sqrt(m_two)
785 deltaf = sqrt(m_e) * xp * exp(-xp * xp)
786
787 case (smear_gaussian)
788 deltaf = exp(-xx * xx) / sqrt(m_pi)
789
790 case (smear_lorentzian)
791 deltaf = m_one / (m_pi * (xx * xx + m_one))
792
793 end select
794
795 end function smear_delta_function
796
797
798 ! ---------------------------------------------------------
799 real(real64) function smear_step_function(this, xx) result(stepf)
800 type(smear_t), intent(in) :: this
801 real(real64), intent(in) :: xx
802
803 real(real64), parameter :: maxarg = 200.0_real64
804 real(real64) :: xp, arg, hd, hp, aa
805 integer :: ii, ni
806
807 push_sub(smear_step_function)
808
809 ! smear_step_function is not defined for SMEAR_FIXED_OCC
810 assert(this%method /= smear_fixed_occ)
811
812 stepf = m_zero
813 select case (this%method)
815 if (xx > m_zero) then
816 stepf = m_one
817 else if (abs(xx) <= m_epsilon) then
818 stepf = this%ef_occ
819 end if
820
821 case (smear_fermi_dirac)
822 if (xx > maxarg) then
823 stepf = m_one
824 else if (xx > -maxarg) then
825 stepf = m_one / (m_one + exp(-xx))
826 end if
827
829 xp = xx - m_one / sqrt(m_two)
830 arg = min(maxarg, xp**2)
831
832 stepf = m_half * erf(xp) + &
833 m_one / sqrt(m_two * m_pi) * exp(-arg) + m_half
834
836 stepf = m_half * erfc(-xx)
837
838 if (this%MP_n > 0) then ! recursion
839 hd = m_zero
840 arg = min(maxarg, xx**2)
841 hp = exp(-arg)
842 ni = 0
843 aa = m_one / sqrt(m_pi)
844 do ii = 1, this%MP_n
845 hd = m_two * xx * hp - m_two * ni * hd
846 ni = ni + 1
847 aa = -aa / (m_four * ii)
848 stepf = stepf - aa * hd
849 hp = m_two * xx * hd - m_two * ni * hp
850 ni = ni + 1
851 end do
852 end if
853
854 case (smear_spline)
855 if (xx <= m_zero) then
856 xp = xx - m_one / sqrt(m_two)
857 stepf = m_half * sqrt(m_e) * exp(-xp * xp)
858 else
859 xp = xx + m_one / sqrt(m_two)
860 stepf = m_one - m_half * sqrt(m_e) * exp(-xp * xp)
861 end if
862
863 case (smear_gaussian)
864 stepf = m_half * erfc(-xx)
865
866 case (smear_lorentzian)
867 stepf = atan(xx) / m_pi + m_half
868
869 end select
870
871 pop_sub(smear_step_function)
872 end function smear_step_function
873
874
875 ! ---------------------------------------------------------
877 real(real64) function smear_entropy_function(this, xx) result(entropyf)
878 type(smear_t), intent(in) :: this
879 real(real64), intent(in) :: xx
880
881 real(real64), parameter :: maxarg = 200.0_real64
882 real(real64) :: xp, arg, hd, hp, hpm1, aa
883 integer :: ii, ni
884
885 push_sub(smear_entropy_function)
886
887 ! smear_entropy_function is not defined for SMEAR_FIXED_OCC
888 assert(this%method /= smear_fixed_occ)
889
890 entropyf = m_zero
891 select case (this%method)
893
895 if (abs(xx) <= 36.0_real64) then
896 xp = m_one / (m_one + exp(-xx))
897 entropyf = xp * log(xp) + (m_one - xp) * log(m_one - xp)
898 end if
899
900 case (smear_cold)
901 xp = xx - m_one / sqrt(m_two)
902 arg = min(maxarg, xp**2)
903
904 entropyf = m_one / sqrt(m_two * m_pi) * xp * exp(-arg)
905
907 arg = min(maxarg, xx**2)
908 entropyf = -m_half * exp(-arg) / sqrt(m_pi)
909
910 if (this%MP_n > 0) then ! recursion
911 hd = m_zero
912 hp = exp(-arg)
913 ni = 0
914 aa = m_one / sqrt(m_pi)
915 do ii = 1, this%MP_n
916 hd = m_two * xx * hp - m_two * ni * hd
917 ni = ni + 1
918 hpm1 = hp
919 hp = m_two * xx * hd - m_two * ni * hp
920 ni = ni + 1
921 aa = -aa / (m_four * ii)
922 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
923 end do
924 end if
925
926 case (smear_spline)
927 xp = abs(xx) + m_one / sqrt(m_two)
928 entropyf = -sqrt(m_e) * (abs(xx) * exp(-xp * xp) / m_two + sqrt(m_pi) / m_four * erfc(xp))
929
930 case (smear_gaussian)
931 entropyf = -exp(-xx * xx) / m_two / sqrt(m_pi)
932
933 case (smear_lorentzian)
934 entropyf = m_ninf
935
936 end select
937
939 end function smear_entropy_function
940
941
942 ! ---------------------------------------------------------
943 logical pure function smear_is_semiconducting(this) result(answer)
944 type(smear_t), intent(in) :: this
945
946 answer = this%method == smear_semiconductor
947
948 end function smear_is_semiconducting
949
950 subroutine smear_write_info(this, namespace, iunit)
951 type(smear_t), intent(in) :: this
952 type(namespace_t), intent(in) :: namespace
953 integer, optional, intent(in) :: iunit
954
955
956 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
957 if (this%photodop) then
958 write(message(1), '(a,f12.6,1x,a)') "Fermi energy (valence ) = ", &
959 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
960 write(message(2), '(a,f12.6,1x,a)') "Fermi energy (conduction) = ", &
961 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
962 call messages_info(2, iunit=iunit, namespace=namespace)
963 else
964 write(message(1), '(a,f12.6,1x,a)') "Fermi energy = ", &
965 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
966 call messages_info(1, iunit=iunit, namespace=namespace)
967 end if
968 end if
969
970 end subroutine smear_write_info
971
973end module smear_oct_m
974
975!! Local Variables:
976!! mode: f90
977!! coding: utf-8
978!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:156
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double atan(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public p_ry
Definition: global.F90:230
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:218
integer function, public kpoints_kweight_denominator(this)
Definition: kpoints.F90:1748
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
recursive integer function, public factorial(n)
Definition: math.F90:289
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public simplex_weights(rdim, esimplex, eF, weights, dos)
Get the weights and DOS contribution of a single simplex.
Definition: simplex.F90:418
integer, parameter, public smear_tetrahedra_opt
Definition: smear.F90:177
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
Definition: smear.F90:781
integer, parameter, public smear_tetrahedra
Definition: smear.F90:177
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
Definition: smear.F90:373
real(real64) function, public smear_entropy_function(this, xx)
This function is defined as .
Definition: smear.F90:973
real(real64) function, public smear_delta_function(this, xx)
Definition: smear.F90:829
integer, parameter, public smear_fermi_dirac
Definition: smear.F90:177
subroutine bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, q_in, e_fermi)
Finds the Fermi energy using the tetrahedron method via bisection.
Definition: smear.F90:585
integer, parameter, public smear_methfessel_paxton
Definition: smear.F90:177
subroutine, public smear_write_info(this, namespace, iunit)
Definition: smear.F90:1046
real(real64) function, public smear_step_function(this, xx)
Definition: smear.F90:895
integer, parameter, public smear_semiconductor
Definition: smear.F90:177
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, kweights, nik, nst)
Definition: smear.F90:668
integer, parameter, public smear_lorentzian
Definition: smear.F90:177
subroutine, public smear_copy(to, from)
Definition: smear.F90:351
integer, parameter, public smear_fixed_occ
Definition: smear.F90:177
integer, parameter, public smear_spline
Definition: smear.F90:177
subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
Definition: smear.F90:512
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
Definition: smear.F90:195
integer, parameter, public smear_cold
Definition: smear.F90:177
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:1039
integer, parameter, public smear_gaussian
Definition: smear.F90:177
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
int true(void)