26 use,
intrinsic :: iso_fortran_env
57 integer,
public :: method
58 real(real64),
public :: dsmear
59 real(real64) :: dsmear_cond
60 real(real64),
public :: e_fermi
61 real(real64) :: e_fermi_cond
63 integer,
public :: el_per_state
64 real(real64),
public :: ef_occ
65 real(real64) :: ef_occ_cond
66 logical,
public :: integral_occs
67 integer,
public :: MP_n
68 integer :: fermi_count
72 type(simplex_t),
pointer :: tetra_simplex => null()
75 logical,
public:: photodop
76 real(real64) :: nephotodop
77 integer :: photodop_bandmin
81 integer,
parameter,
public :: &
82 SMEAR_SEMICONDUCTOR = 1, &
93 real(real64),
parameter :: TOL_SMEAR = 1e-6_real64
98 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
99 type(smear_t),
intent(out) :: this
100 type(namespace_t),
intent(in) :: namespace
101 integer,
intent(in) :: ispin
102 logical,
intent(in) :: fixed_occ
103 logical,
intent(in) :: integral_occs
104 type(kpoints_t),
intent(in) :: kpoints
108 this%integral_occs = integral_occs
146 call parse_variable(namespace,
'SmearingFunction', smear_semiconductor, this%method)
160 this%dsmear = 1e-14_real64
174 this%dsmear_cond = this%dsmear
175 if (this%method /= smear_semiconductor .and. this%method /=
smear_fixed_occ)
then
192 this%photodop=.false.
203 call parse_variable(namespace,
'PhotoDopingBand', 1, this%photodop_bandmin)
209 this%el_per_state = 2
211 this%el_per_state = 1
220 if (this%method == smear_semiconductor)
then
223 if (this%nik_factor == 0)
then
224 message(1) =
"k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
228 if (.not.
associated(kpoints%reduced%simplex))
then
229 message(1) =
"Tetrahedron smearing requires a simplex mesh for the k-point grid."
233 this%tetra_simplex => kpoints%reduced%simplex
254 type(
smear_t),
intent(out) :: to
255 type(
smear_t),
intent(in) :: from
259 to%method = from%method
260 to%dsmear = from%dsmear
261 to%e_fermi = from%e_fermi
262 to%el_per_state = from%el_per_state
263 to%ef_occ = from%ef_occ
265 to%fermi_count = from%fermi_count
266 to%nik_factor = from%nik_factor
267 to%tetra_simplex => from%tetra_simplex
275 qtot, nik, nst, kweights)
276 type(
smear_t),
intent(inout) :: this
278 real(real64),
intent(in) :: eigenvalues(:,:), occupations(:,:)
279 real(real64),
intent(in) :: qtot, kweights(:)
280 integer,
intent(in) :: nik, nst
282 real(real64),
parameter :: tol = 1.0e-10_real64
283 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
284 real(real64) :: sumq_frac
286 real(real64),
allocatable :: eigenval_list(:)
287 integer,
allocatable :: k_list(:), reorder(:)
288 integer :: fermi_count_up, fermi_count_down
292 maxq = this%el_per_state * nst * this%nspins
293 if (maxq - qtot <= -tol)
then
294 message(1) =
'Not enough states'
295 write(
message(2),
'(6x,a,f12.6,a,i10)')
'(total charge = ', qtot, &
296 ' max charge = ', maxq
302 ist_cycle:
do ist = nst, 1, -1
303 if (any(occupations(ist, :) >
m_min_occ))
then
304 this%e_fermi = eigenvalues(ist, 1)
305 this%ef_occ = occupations(ist, 1) / this%el_per_state
307 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) >
m_min_occ)
then
308 this%e_fermi = eigenvalues(ist, ik)
309 this%ef_occ = occupations(ist, ik) / this%el_per_state
316 else if (this%method == smear_semiconductor)
then
318 safe_allocate(eigenval_list(1:nst * nik))
319 safe_allocate( k_list(1:nst * nik))
320 safe_allocate( reorder(1:nst * nik))
325 eigenval_list(iter) = eigenvalues(ist, ik)
332 call sort(eigenval_list, reorder)
334 sumq_int =
floor(qtot) * this%nik_factor
335 sumq_frac = qtot * this%nik_factor - sumq_int
336 if ( sumq_frac + tol > this%nik_factor )
then
337 sumq_int = sumq_int + this%nik_factor
338 sumq_frac = sumq_frac - this%nik_factor
340 if (sumq_frac < tol) sumq_frac =
m_zero
342 do iter = 1, nst * nik
343 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor +
m_half)
344 if (.not. weight > 0) cycle
345 this%e_fermi = eigenval_list(iter)
346 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
348 if (sumq_int - weight * this%el_per_state <= 0)
then
356 if (iter - fermi_count_down < 1)
exit
357 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear)
exit
358 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor +
m_half)
360 sumq_int = sumq_int + weight * this%el_per_state
361 sum_weight = sum_weight + weight
363 fermi_count_down = fermi_count_down + 1
366 if (iter + fermi_count_up > nst*nik)
exit
367 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear)
exit
368 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor +
m_half)
370 sum_weight = sum_weight + weight
372 fermi_count_up = fermi_count_up + 1
374 this%fermi_count = fermi_count_up + fermi_count_down - 1
375 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
379 sumq_int = sumq_int - weight * this%el_per_state
383 safe_deallocate_a(eigenval_list)
384 safe_deallocate_a(k_list)
385 safe_deallocate_a(reorder)
388 if (this%photodop)
then
393 if (this%photodop)
then
396 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
397 this%e_fermi, this%ef_occ)
401 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
402 this%e_fermi_cond, this%ef_occ_cond)
405 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
414 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
415 type(
smear_t),
intent(inout) :: this
417 real(real64),
intent(in) :: dsmear_in
418 real(real64),
intent(in) :: tol
419 integer,
intent(in) :: nik
420 real(real64),
intent(in) :: eigenvalues(:,:)
421 real(real64),
intent(in) :: kweights(:), q_in
422 integer,
intent(in) :: start_band, end_band
423 real(real64),
intent(out) :: e_fermi, ef_occ
425 integer,
parameter :: nitmax = 200
426 integer :: ist, ik, iter
427 real(real64) :: drange, dsmear, emin, emax, xx, sumq
432 dsmear = max(1e-14_real64, dsmear_in)
433 drange = dsmear *
sqrt(-
log(tol * 0.01_real64))
435 emin = minval(eigenvalues) - drange
436 emax = maxval(eigenvalues) + drange
439 e_fermi =
m_half * (emin + emax)
443 do ist = start_band, end_band
444 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
445 sumq = sumq + kweights(ik) * this%el_per_state * &
450 conv = (abs(sumq - q_in) <= tol)
453 if (sumq <= q_in) emin = e_fermi
454 if (sumq >= q_in) emax = e_fermi
460 message(1) =
'Fermi: did not converge.'
488 type(
smear_t),
intent(inout) :: this
490 real(real64),
intent(in) :: tol
491 integer,
intent(in) :: nik, nst
492 real(real64),
intent(in) :: eigenvalues(:,:)
493 real(real64),
intent(in) :: q_in
494 real(real64),
intent(out) :: e_fermi
496 integer,
parameter :: nitmax = 200
498 real(real64) :: emin, emax, sumq
501 integer :: ist, ii, ll, is, ik, ns, nik_red, nvertices
502 real(real64) :: e_simplex(20), w_simplex(4), dos_simplex(4)
506 assert(
associated(this%tetra_simplex))
507 assert(this%tetra_simplex%n_simplices > 0)
508 assert(this%tetra_simplex%sdim > 0)
512 assert(mod(nik, ns) == 0)
515 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
517 nvertices = this%tetra_simplex%rdim + 1
518 assert(nvertices > 0)
519 assert(nvertices <= this%tetra_simplex%sdim)
521 emin = minval(eigenvalues)
522 emax = maxval(eigenvalues)
523 if (.not. (emin < emax))
then
524 write(
message(1),
'(a)')
'Tetra bisection: invalid energy bracket from eigenvalues.'
529 e_fermi =
m_half * (emin + emax)
533 do ii = 1, this%tetra_simplex%n_simplices
535 do ll = 1, this%tetra_simplex%sdim
536 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
537 e_simplex(ll) = eigenvalues(ist, ik + is)
541 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), e_fermi, &
542 w_simplex(1:nvertices), dos_simplex(1:nvertices))
545 sumq = sumq + this%el_per_state * w_simplex(ll) / this%tetra_simplex%n_points
552 conv = (abs(sumq - q_in) <= tol)
555 if (sumq <= q_in) emin = e_fermi
556 if (sumq >= q_in) emax = e_fermi
560 message(1) =
'Fermi (tetrahedra): did not converge.'
569 type(
smear_t),
intent(in) :: this
570 real(real64),
intent(in) :: eigenvalues(:,:)
571 real(real64),
intent(inout) :: occupations(:,:)
572 real(real64),
intent(in) :: kweights(:)
573 integer,
intent(in) :: nik, nst
575 integer :: ik, ist, ifermi
576 real(real64) :: dsmear, xx, dsmear_cond
582 else if (this%method == smear_semiconductor)
then
583 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
588 xx = eigenvalues(ist, ik) - this%e_fermi
589 if (xx < -tol_smear)
then
590 occupations(ist, ik) = this%el_per_state
591 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count)
then
592 occupations(ist, ik) = this%ef_occ * this%el_per_state
595 occupations(ist, ik) =
m_zero
601 if (this%photodop)
then
605 integer :: ii, ll, is, ns, nik_red, nvertices
606 real(real64) :: e_simplex(20), w_simplex(4), dos_simplex(4)
608 assert(
associated(this%tetra_simplex))
609 assert(this%tetra_simplex%n_simplices > 0)
610 assert(this%tetra_simplex%sdim > 0)
614 assert(mod(nik, ns) == 0)
617 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
619 nvertices = this%tetra_simplex%rdim + 1
620 assert(nvertices > 0)
621 assert(nvertices <= this%tetra_simplex%sdim)
623 occupations(:, :) =
m_zero
626 do ii = 1, this%tetra_simplex%n_simplices
628 do ll = 1, this%tetra_simplex%sdim
629 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
630 e_simplex(ll) = eigenvalues(ist, ik + is)
634 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), this%e_fermi, &
635 w_simplex(1:nvertices), dos_simplex(1:nvertices))
638 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
640 occupations(ist, ik + is) = occupations(ist, ik + is) + this%el_per_state * &
641 w_simplex(ll) / (this%tetra_simplex%n_points * kweights(ik + is))
649 dsmear = max(1e-14_real64, this%dsmear)
650 if (this%photodop)
then
651 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
654 if (ist < this%photodop_bandmin)
then
656 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
659 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
667 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
679 nik, nst, kweights, occ) result(entropy)
680 type(
smear_t),
intent(inout) :: this
681 real(real64),
intent(in) :: eigenvalues(:,:)
682 real(real64),
intent(in) :: kweights(:)
683 integer,
intent(in) :: nik, nst
684 real(real64),
intent(in) :: occ(:, :)
687 real(real64) :: dsmear, xx, term, ff
699 ff = occ(ist, ik) / this%el_per_state
700 if (ff > m_zero .and. ff < m_one)
then
701 term = ff *
log(ff) + (1 - ff) *
log(1 - ff)
705 entropy = entropy - kweights(ik) * this%el_per_state * term
709 dsmear = max(1e-14_real64, this%dsmear)
713 if (eigenvalues(ist, ik) < huge(m_one))
then
714 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
715 entropy = entropy - kweights(ik) * this%el_per_state * &
728 type(
smear_t),
intent(in) :: this
729 real(real64),
intent(in) :: xx
731 real(real64),
parameter :: maxarg = 200.0_real64
732 real(real64) :: xp, arg, hd, hp, aa
741 select case (this%method)
743 if (abs(xx) <= m_epsilon)
then
748 if (abs(xx) <= 36.0_real64)
then
749 deltaf = m_one / (m_two +
exp(-xx) +
exp(xx))
753 xp = xx - m_one /
sqrt(m_two)
754 arg = min(maxarg, xp**2)
756 deltaf =
exp(-arg) /
sqrt(m_pi) * (m_two -
sqrt(m_two) * xx)
759 arg = min(maxarg, xx**2)
760 deltaf =
exp(-arg) /
sqrt(m_pi)
762 if (this%MP_n > 0)
then
766 aa = m_one /
sqrt(m_pi)
768 hd = m_two * xx * hp - m_two * ni * hd
770 aa = -aa / (m_four * ii)
771 hp = m_two * xx * hd - m_two * ni * hp
773 deltaf = deltaf + aa * hp
778 xp = abs(xx) + m_one /
sqrt(m_two)
779 deltaf =
sqrt(m_e) * xp *
exp(-xp * xp)
782 deltaf =
exp(-xx * xx) /
sqrt(m_pi)
785 deltaf = m_one / (m_pi * (xx * xx + m_one))
794 type(
smear_t),
intent(in) :: this
795 real(real64),
intent(in) :: xx
797 real(real64),
parameter :: maxarg = 200.0_real64
798 real(real64) :: xp, arg, hd, hp, aa
807 select case (this%method)
809 if (xx > m_zero)
then
811 else if (abs(xx) <= m_epsilon)
then
816 if (xx > maxarg)
then
818 else if (xx > -maxarg)
then
819 stepf = m_one / (m_one +
exp(-xx))
823 xp = xx - m_one /
sqrt(m_two)
824 arg = min(maxarg, xp**2)
826 stepf = m_half * erf(xp) + &
827 m_one /
sqrt(m_two * m_pi) *
exp(-arg) + m_half
830 stepf = m_half * erfc(-xx)
832 if (this%MP_n > 0)
then
834 arg = min(maxarg, xx**2)
837 aa = m_one /
sqrt(m_pi)
839 hd = m_two * xx * hp - m_two * ni * hd
841 aa = -aa / (m_four * ii)
842 stepf = stepf - aa * hd
843 hp = m_two * xx * hd - m_two * ni * hp
849 if (xx <= m_zero)
then
850 xp = xx - m_one /
sqrt(m_two)
851 stepf = m_half *
sqrt(m_e) *
exp(-xp * xp)
853 xp = xx + m_one /
sqrt(m_two)
854 stepf = m_one - m_half *
sqrt(m_e) *
exp(-xp * xp)
858 stepf = m_half * erfc(-xx)
861 stepf =
atan(xx) / m_pi + m_half
872 type(
smear_t),
intent(in) :: this
873 real(real64),
intent(in) :: xx
875 real(real64),
parameter :: maxarg = 200.0_real64
876 real(real64) :: xp, arg, hd, hp, hpm1, aa
885 select case (this%method)
889 if (abs(xx) <= 36.0_real64)
then
890 xp = m_one / (m_one +
exp(-xx))
891 entropyf = xp *
log(xp) + (m_one - xp) *
log(m_one - xp)
895 xp = xx - m_one /
sqrt(m_two)
896 arg = min(maxarg, xp**2)
898 entropyf = m_one /
sqrt(m_two * m_pi) * xp *
exp(-arg)
901 arg = min(maxarg, xx**2)
902 entropyf = -m_half *
exp(-arg) /
sqrt(m_pi)
904 if (this%MP_n > 0)
then
908 aa = m_one /
sqrt(m_pi)
910 hd = m_two * xx * hp - m_two * ni * hd
913 hp = m_two * xx * hd - m_two * ni * hp
915 aa = -aa / (m_four * ii)
916 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
921 xp = abs(xx) + m_one /
sqrt(m_two)
922 entropyf = -
sqrt(m_e) * (abs(xx) *
exp(-xp * xp) / m_two +
sqrt(m_pi) / m_four * erfc(xp))
925 entropyf = -
exp(-xx * xx) / m_two /
sqrt(m_pi)
938 type(
smear_t),
intent(in) :: this
945 type(
smear_t),
intent(in) :: this
946 type(namespace_t),
intent(in) :: namespace
947 integer,
optional,
intent(in) :: iunit
951 if (this%photodop)
then
952 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy (valence ) = ", &
953 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
954 write(message(2),
'(a,f12.6,1x,a)')
"Fermi energy (conduction) = ", &
955 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
956 call messages_info(2, iunit=iunit, namespace=namespace)
958 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy = ", &
959 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
960 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)
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