26 use,
intrinsic :: iso_fortran_env
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
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
73 type(simplex_t),
pointer :: tetra_simplex => null()
76 logical,
public:: photodop
77 real(real64) :: nephotodop
78 integer :: photodop_bandmin
82 integer,
parameter,
public :: &
83 SMEAR_SEMICONDUCTOR = 1, &
94 real(real64),
parameter :: TOL_SMEAR = 1e-6_real64
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
109 this%integral_occs = integral_occs
148 call parse_variable(namespace,
'SmearingFunction', smear_semiconductor, this%method)
162 this%dsmear = 1e-14_real64
176 this%dsmear_cond = this%dsmear
194 this%photodop=.false.
205 call parse_variable(namespace,
'PhotoDopingBand', 1, this%photodop_bandmin)
211 this%el_per_state = 2
213 this%el_per_state = 1
222 if (this%method == smear_semiconductor)
then
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."
230 if (.not.
associated(kpoints%simplex))
then
231 message(1) =
"Tetrahedron smearing requires a simplex mesh for the k-point grid."
235 this%tetra_simplex => kpoints%simplex
256 type(
smear_t),
intent(out) :: to
257 type(
smear_t),
intent(in) :: from
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
267 to%fermi_count = from%fermi_count
268 to%nik_factor = from%nik_factor
269 to%tetra_simplex => from%tetra_simplex
277 qtot, nik, nst, kweights)
278 type(
smear_t),
intent(inout) :: this
280 real(real64),
intent(in) :: eigenvalues(:,:), occupations(:,:)
281 real(real64),
intent(in) :: qtot, kweights(:)
282 integer,
intent(in) :: nik, nst
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
288 real(real64),
allocatable :: eigenval_list(:)
289 integer,
allocatable :: k_list(:), reorder(:)
290 integer :: fermi_count_up, fermi_count_down
294 maxq = this%el_per_state * nst * this%nspins
295 if (maxq - qtot <= -tol)
then
296 message(1) =
'Not enough states'
297 write(
message(2),
'(6x,a,f12.6,a,i10)')
'(total charge = ', qtot, &
298 ' max charge = ', maxq
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
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
318 else if (this%method == smear_semiconductor)
then
320 safe_allocate(eigenval_list(1:nst * nik))
321 safe_allocate( k_list(1:nst * nik))
322 safe_allocate( reorder(1:nst * nik))
327 eigenval_list(iter) = eigenvalues(ist, ik)
334 call sort(eigenval_list, reorder)
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
342 if (sumq_frac < tol) sumq_frac =
m_zero
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)
350 if (sumq_int - weight * this%el_per_state <= 0)
then
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)
362 sumq_int = sumq_int + weight * this%el_per_state
363 sum_weight = sum_weight + weight
365 fermi_count_down = fermi_count_down + 1
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)
372 sum_weight = sum_weight + weight
374 fermi_count_up = fermi_count_up + 1
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)
381 sumq_int = sumq_int - weight * this%el_per_state
385 safe_deallocate_a(eigenval_list)
386 safe_deallocate_a(k_list)
387 safe_deallocate_a(reorder)
390 if (this%photodop)
then
395 if (this%photodop)
then
398 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
399 this%e_fermi, this%ef_occ)
403 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
404 this%e_fermi_cond, this%ef_occ_cond)
407 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
416 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
417 type(
smear_t),
intent(inout) :: this
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
427 integer,
parameter :: nitmax = 200
428 integer :: ist, ik, iter
429 real(real64) :: drange, dsmear, emin, emax, xx, sumq
434 dsmear = max(1e-14_real64, dsmear_in)
435 drange = dsmear *
sqrt(-
log(tol * 0.01_real64))
437 emin = minval(eigenvalues) - drange
438 emax = maxval(eigenvalues) + drange
441 e_fermi =
m_half * (emin + emax)
445 do ist = start_band, end_band
446 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
447 sumq = sumq + kweights(ik) * this%el_per_state * &
452 conv = (abs(sumq - q_in) <= tol)
455 if (sumq <= q_in) emin = e_fermi
456 if (sumq >= q_in) emax = e_fermi
462 message(1) =
'Fermi: did not converge.'
490 type(
smear_t),
intent(inout) :: this
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
498 integer,
parameter :: nitmax = 200
500 real(real64) :: emin, emax, sumq
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
508 assert(
associated(this%tetra_simplex))
509 assert(this%tetra_simplex%n_simplices > 0)
510 assert(this%tetra_simplex%sdim > 0)
514 assert(mod(nik, ns) == 0)
517 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
519 nvertices = this%tetra_simplex%rdim + 1
520 assert(nvertices > 0)
521 assert(nvertices <= this%tetra_simplex%sdim)
523 simplex_unit_volume =
factorial(this%tetra_simplex%rdim)
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.'
533 e_fermi =
m_half * (emin + emax)
537 do ii = 1, this%tetra_simplex%n_simplices
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)
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)
549 sumq = sumq + this%el_per_state * simplex_unit_volume * w_simplex(ll) / this%tetra_simplex%n_simplices
556 conv = (abs(sumq - q_in) <= tol)
559 if (sumq <= q_in) emin = e_fermi
560 if (sumq >= q_in) emax = e_fermi
564 message(1) =
'Fermi (tetrahedra): did not converge.'
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
579 integer :: ik, ist, ifermi
580 real(real64) :: dsmear, xx, dsmear_cond
586 else if (this%method == smear_semiconductor)
then
587 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*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
599 occupations(ist, ik) =
m_zero
605 if (this%photodop)
then
609 integer :: ii, ll, is, ns, nik_red, nvertices
610 real(real64) :: e_simplex(20), w_simplex(20), dos_simplex, simplex_unit_volume
612 assert(
associated(this%tetra_simplex))
613 assert(this%tetra_simplex%n_simplices > 0)
614 assert(this%tetra_simplex%sdim > 0)
618 assert(mod(nik, ns) == 0)
621 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
623 nvertices = this%tetra_simplex%rdim + 1
624 assert(nvertices > 0)
625 assert(nvertices <= this%tetra_simplex%sdim)
627 simplex_unit_volume =
factorial(this%tetra_simplex%rdim)
629 occupations(:, :) =
m_zero
632 do ii = 1, this%tetra_simplex%n_simplices
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)
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)
644 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
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))
655 dsmear = max(1e-14_real64, this%dsmear)
656 if (this%photodop)
then
657 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
660 if (ist < this%photodop_bandmin)
then
662 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
665 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
673 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
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(:, :)
693 real(real64) :: dsmear, xx, term, ff
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)
711 entropy = entropy - kweights(ik) * this%el_per_state * term
715 dsmear = max(1e-14_real64, this%dsmear)
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 * &
734 type(
smear_t),
intent(in) :: this
735 real(real64),
intent(in) :: xx
737 real(real64),
parameter :: maxarg = 200.0_real64
738 real(real64) :: xp, arg, hd, hp, aa
747 select case (this%method)
749 if (abs(xx) <= m_epsilon)
then
754 if (abs(xx) <= 36.0_real64)
then
755 deltaf = m_one / (m_two +
exp(-xx) +
exp(xx))
759 xp = xx - m_one /
sqrt(m_two)
760 arg = min(maxarg, xp**2)
762 deltaf =
exp(-arg) /
sqrt(m_pi) * (m_two -
sqrt(m_two) * xx)
765 arg = min(maxarg, xx**2)
766 deltaf =
exp(-arg) /
sqrt(m_pi)
768 if (this%MP_n > 0)
then
772 aa = m_one /
sqrt(m_pi)
774 hd = m_two * xx * hp - m_two * ni * hd
776 aa = -aa / (m_four * ii)
777 hp = m_two * xx * hd - m_two * ni * hp
779 deltaf = deltaf + aa * hp
784 xp = abs(xx) + m_one /
sqrt(m_two)
785 deltaf =
sqrt(m_e) * xp *
exp(-xp * xp)
788 deltaf =
exp(-xx * xx) /
sqrt(m_pi)
791 deltaf = m_one / (m_pi * (xx * xx + m_one))
800 type(
smear_t),
intent(in) :: this
801 real(real64),
intent(in) :: xx
803 real(real64),
parameter :: maxarg = 200.0_real64
804 real(real64) :: xp, arg, hd, hp, aa
813 select case (this%method)
815 if (xx > m_zero)
then
817 else if (abs(xx) <= m_epsilon)
then
822 if (xx > maxarg)
then
824 else if (xx > -maxarg)
then
825 stepf = m_one / (m_one +
exp(-xx))
829 xp = xx - m_one /
sqrt(m_two)
830 arg = min(maxarg, xp**2)
832 stepf = m_half * erf(xp) + &
833 m_one /
sqrt(m_two * m_pi) *
exp(-arg) + m_half
836 stepf = m_half * erfc(-xx)
838 if (this%MP_n > 0)
then
840 arg = min(maxarg, xx**2)
843 aa = m_one /
sqrt(m_pi)
845 hd = m_two * xx * hp - m_two * ni * hd
847 aa = -aa / (m_four * ii)
848 stepf = stepf - aa * hd
849 hp = m_two * xx * hd - m_two * ni * hp
855 if (xx <= m_zero)
then
856 xp = xx - m_one /
sqrt(m_two)
857 stepf = m_half *
sqrt(m_e) *
exp(-xp * xp)
859 xp = xx + m_one /
sqrt(m_two)
860 stepf = m_one - m_half *
sqrt(m_e) *
exp(-xp * xp)
864 stepf = m_half * erfc(-xx)
867 stepf =
atan(xx) / m_pi + m_half
878 type(
smear_t),
intent(in) :: this
879 real(real64),
intent(in) :: xx
881 real(real64),
parameter :: maxarg = 200.0_real64
882 real(real64) :: xp, arg, hd, hp, hpm1, aa
891 select case (this%method)
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)
901 xp = xx - m_one /
sqrt(m_two)
902 arg = min(maxarg, xp**2)
904 entropyf = m_one /
sqrt(m_two * m_pi) * xp *
exp(-arg)
907 arg = min(maxarg, xx**2)
908 entropyf = -m_half *
exp(-arg) /
sqrt(m_pi)
910 if (this%MP_n > 0)
then
914 aa = m_one /
sqrt(m_pi)
916 hd = m_two * xx * hp - m_two * ni * hd
919 hp = m_two * xx * hd - m_two * ni * hp
921 aa = -aa / (m_four * ii)
922 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
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))
931 entropyf = -
exp(-xx * xx) / m_two /
sqrt(m_pi)
944 type(
smear_t),
intent(in) :: this
951 type(
smear_t),
intent(in) :: this
952 type(namespace_t),
intent(in) :: namespace
953 integer,
optional,
intent(in) :: iunit
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)
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)
This is the common interface to a sorting routine. It performs the shell algorithm,...
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
real(real64), parameter, public m_zero
real(real64), parameter, public p_ry
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
integer function, public kpoints_kweight_denominator(this)
This module is intended to contain "only mathematical" functions and procedures.
recursive integer function, public factorial(n)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public simplex_weights(rdim, esimplex, eF, weights, dos)
Get the weights and DOS contribution of a single simplex.
integer, parameter, public smear_tetrahedra_opt
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
integer, parameter, public smear_tetrahedra
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
real(real64) function, public smear_entropy_function(this, xx)
This function is defined as .
real(real64) function, public smear_delta_function(this, xx)
integer, parameter, public smear_fermi_dirac
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.
integer, parameter, public smear_methfessel_paxton
subroutine, public smear_write_info(this, namespace, iunit)
real(real64) function, public smear_step_function(this, xx)
integer, parameter, public smear_semiconductor
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, kweights, nik, nst)
integer, parameter, public smear_lorentzian
subroutine, public smear_copy(to, from)
integer, parameter, public smear_fixed_occ
integer, parameter, public smear_spline
subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
integer, parameter, public smear_cold
logical pure function, public smear_is_semiconducting(this)
integer, parameter, public smear_gaussian
This module is intended to contain "only mathematical" functions and procedures.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
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