26 use,
intrinsic :: iso_fortran_env
56 integer,
public :: method
57 real(real64),
public :: dsmear
58 real(real64) :: dsmear_cond
59 real(real64),
public :: e_fermi
60 real(real64) :: e_fermi_cond
62 integer,
public :: el_per_state
63 real(real64),
public :: ef_occ
64 real(real64) :: ef_occ_cond
65 logical,
public :: integral_occs
67 integer :: fermi_count
72 logical,
public:: photodop
73 real(real64) :: nephotodop
74 integer :: photodop_bandmin
78 integer,
parameter,
public :: &
79 SMEAR_SEMICONDUCTOR = 1, &
86 real(real64),
parameter :: TOL_SMEAR = 1e-6_real64
91 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
92 type(smear_t),
intent(out) :: this
93 type(namespace_t),
intent(in) :: namespace
94 integer,
intent(in) :: ispin
95 logical,
intent(in) :: fixed_occ
96 logical,
intent(in) :: integral_occs
97 type(kpoints_t),
intent(in) :: kpoints
101 this%integral_occs = integral_occs
129 call parse_variable(namespace,
'SmearingFunction', smear_semiconductor, this%method)
143 this%dsmear = 1e-14_real64
144 if (this%method /= smear_semiconductor .and. this%method /=
smear_fixed_occ)
then
157 this%dsmear_cond = this%dsmear
175 this%photodop=.false.
186 call parse_variable(namespace,
'PhotoDopingBand', 1, this%photodop_bandmin)
192 this%el_per_state = 2
194 this%el_per_state = 1
203 if (this%method == smear_semiconductor)
then
206 if (this%nik_factor == 0)
then
207 message(1) =
"k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
230 type(
smear_t),
intent(out) :: to
231 type(
smear_t),
intent(in) :: from
235 to%method = from%method
236 to%dsmear = from%dsmear
237 to%e_fermi = from%e_fermi
238 to%el_per_state = from%el_per_state
239 to%ef_occ = from%ef_occ
241 to%fermi_count = from%fermi_count
242 to%nik_factor = from%nik_factor
244 to%dsmear_cond = from%dsmear_cond
245 to%e_fermi_cond = from%e_fermi_cond
246 to%ef_occ_cond = from%ef_occ_cond
247 to%integral_occs = from%integral_occs
248 to%nspins = from%nspins
249 to%photodop = from%photodop
250 to%nephotodop = from%nephotodop
251 to%photodop_bandmin = from%photodop_bandmin
259 qtot, nik, nst, kweights)
260 type(
smear_t),
intent(inout) :: this
262 real(real64),
intent(in) :: eigenvalues(:,:), occupations(:,:)
263 real(real64),
intent(in) :: qtot, kweights(:)
264 integer,
intent(in) :: nik, nst
266 real(real64),
parameter :: tol = 1.0e-10_real64
267 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
268 real(real64) :: sumq_frac
270 real(real64),
allocatable :: eigenval_list(:)
271 integer,
allocatable :: k_list(:), reorder(:)
272 integer :: fermi_count_up, fermi_count_down
276 maxq = this%el_per_state * nst * this%nspins
277 if (maxq - qtot <= -tol)
then
278 message(1) =
'Not enough states'
279 write(
message(2),
'(6x,a,f12.6,a,i10)')
'(total charge = ', qtot, &
280 ' max charge = ', maxq
286 ist_cycle:
do ist = nst, 1, -1
287 if (any(occupations(ist, :) >
m_min_occ))
then
288 this%e_fermi = eigenvalues(ist, 1)
289 this%ef_occ = occupations(ist, 1) / this%el_per_state
291 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) >
m_min_occ)
then
292 this%e_fermi = eigenvalues(ist, ik)
293 this%ef_occ = occupations(ist, ik) / this%el_per_state
300 else if (this%method == smear_semiconductor)
then
302 safe_allocate(eigenval_list(1:nst * nik))
303 safe_allocate( k_list(1:nst * nik))
304 safe_allocate( reorder(1:nst * nik))
309 eigenval_list(iter) = eigenvalues(ist, ik)
316 call sort(eigenval_list, reorder)
318 sumq_int =
floor(qtot) * this%nik_factor
319 sumq_frac = qtot * this%nik_factor - sumq_int
320 if ( sumq_frac + tol > this%nik_factor )
then
321 sumq_int = sumq_int + this%nik_factor
322 sumq_frac = sumq_frac - this%nik_factor
324 if (sumq_frac < tol) sumq_frac =
m_zero
326 do iter = 1, nst * nik
327 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor +
m_half)
328 if (weight <= 0) cycle
329 this%e_fermi = eigenval_list(iter)
330 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
332 if (sumq_int - weight * this%el_per_state <= 0)
then
340 if (iter - fermi_count_down < 1)
exit
341 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear)
exit
342 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor +
m_half)
344 sumq_int = sumq_int + weight * this%el_per_state
345 sum_weight = sum_weight + weight
347 fermi_count_down = fermi_count_down + 1
350 if (iter + fermi_count_up > nst*nik)
exit
351 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear)
exit
352 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor +
m_half)
354 sum_weight = sum_weight + weight
356 fermi_count_up = fermi_count_up + 1
358 this%fermi_count = fermi_count_up + fermi_count_down - 1
359 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
363 sumq_int = sumq_int - weight * this%el_per_state
367 safe_deallocate_a(eigenval_list)
368 safe_deallocate_a(k_list)
369 safe_deallocate_a(reorder)
372 if (this%photodop)
then
375 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
376 this%e_fermi, this%ef_occ)
380 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
381 this%e_fermi_cond, this%ef_occ_cond)
384 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
393 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
394 type(
smear_t),
intent(inout) :: this
396 real(real64),
intent(in) :: dsmear_in
397 real(real64),
intent(in) :: tol
398 integer,
intent(in) :: nik
399 real(real64),
intent(in) :: eigenvalues(:,:)
400 real(real64),
intent(in) :: kweights(:), q_in
401 integer,
intent(in) :: start_band, end_band
402 real(real64),
intent(out) :: e_fermi, ef_occ
404 integer,
parameter :: nitmax = 200
405 integer :: ist, ik, iter
406 real(real64) :: drange, dsmear, emin, emax, xx, sumq
411 dsmear = max(1e-14_real64, dsmear_in)
412 drange = dsmear *
sqrt(-
log(tol * 0.01_real64))
414 emin = minval(eigenvalues) - drange
415 emax = maxval(eigenvalues) + drange
418 e_fermi =
m_half * (emin + emax)
422 do ist = start_band, end_band
423 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
424 sumq = sumq + kweights(ik) * this%el_per_state * &
429 conv = (abs(sumq - q_in) <= tol)
432 if (sumq <= q_in) emin = e_fermi
433 if (sumq >= q_in) emax = e_fermi
439 message(1) =
'Fermi: did not converge.'
448 type(
smear_t),
intent(in) :: this
449 real(real64),
intent(in) :: eigenvalues(:,:)
450 real(real64),
intent(inout) :: occupations(:,:)
451 integer,
intent(in) :: nik, nst
453 integer :: ik, ist, ifermi
454 real(real64) :: dsmear, xx, dsmear_cond
460 else if (this%method == smear_semiconductor)
then
461 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
466 xx = eigenvalues(ist, ik) - this%e_fermi
467 if (xx < -tol_smear)
then
468 occupations(ist, ik) = this%el_per_state
469 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count)
then
470 occupations(ist, ik) = this%ef_occ * this%el_per_state
473 occupations(ist, ik) =
m_zero
480 dsmear = max(1e-14_real64, this%dsmear)
481 if (this%photodop)
then
482 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
485 if (ist < this%photodop_bandmin)
then
487 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
490 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
498 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
510 real(real64) function smear_calc_entropy(this, eigenvalues, &
511 nik, nst, kweights, occ) result(entropy)
512 type(
smear_t),
intent(in) :: this
513 real(real64),
intent(in) :: eigenvalues(:,:)
514 real(real64),
intent(in) :: kweights(:)
515 integer,
intent(in) :: nik, nst
516 real(real64),
intent(in) :: occ(:, :)
519 real(real64) :: dsmear, xx, term, ff
521 push_sub(smear_calc_entropy)
525 if (this%method ==
smear_fixed_occ .or. this%method == smear_semiconductor)
then
531 ff = occ(ist, ik) / this%el_per_state
533 term = ff *
log(ff) + (1 - ff) *
log(1 - ff)
537 entropy = entropy - kweights(ik) * this%el_per_state * term
541 dsmear = max(1e-14_real64, this%dsmear)
545 if (eigenvalues(ist, ik) < huge(
m_one))
then
546 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
547 entropy = entropy - kweights(ik) * this%el_per_state * &
554 pop_sub(smear_calc_entropy)
560 type(
smear_t),
intent(in) :: this
561 real(real64),
intent(in) :: xx
563 real(real64),
parameter :: maxarg = 200.0_real64
564 real(real64) :: xp, arg, hd, hp, aa
573 select case (this%method)
575 if (abs(xx) <= m_epsilon)
then
580 if (abs(xx) <= 36.0_real64)
then
581 deltaf = m_one / (m_two +
exp(-xx) +
exp(xx))
585 xp = xx - m_one /
sqrt(m_two)
586 arg = min(maxarg, xp**2)
588 deltaf =
exp(-arg) /
sqrt(m_pi) * (m_two -
sqrt(m_two) * xx)
591 arg = min(maxarg, xx**2)
592 deltaf =
exp(-arg) /
sqrt(m_pi)
594 if (this%MP_n > 0)
then
598 aa = m_one /
sqrt(m_pi)
600 hd = m_two * xx * hp - m_two * ni * hd
602 aa = -aa / (m_four * ii)
603 hp = m_two * xx * hd - m_two * ni * hp
605 deltaf = deltaf + aa * hp
610 xp = abs(xx) + m_one /
sqrt(m_two)
611 deltaf =
sqrt(m_e) * xp *
exp(-xp * xp)
620 type(
smear_t),
intent(in) :: this
621 real(real64),
intent(in) :: xx
623 real(real64),
parameter :: maxarg = 200.0_real64
624 real(real64) :: xp, arg, hd, hp, aa
633 select case (this%method)
635 if (xx > m_zero)
then
637 else if (abs(xx) <= m_epsilon)
then
642 if (xx > maxarg)
then
644 else if (xx > -maxarg)
then
645 stepf = m_one / (m_one +
exp(-xx))
649 xp = xx - m_one /
sqrt(m_two)
650 arg = min(maxarg, xp**2)
652 stepf = m_half * loct_erf(xp) + &
653 m_one /
sqrt(m_two * m_pi) *
exp(-arg) + m_half
656 stepf = m_half * loct_erfc(-xx)
658 if (this%MP_n > 0)
then
660 arg = min(maxarg, xx**2)
663 aa = m_one /
sqrt(m_pi)
665 hd = m_two * xx * hp - m_two * ni * hd
667 aa = -aa / (m_four * ii)
668 stepf = stepf - aa * hd
669 hp = m_two * xx * hd - m_two * ni * hp
675 if (xx <= m_zero)
then
676 xp = xx - m_one /
sqrt(m_two)
677 stepf = m_half *
sqrt(m_e) *
exp(-xp * xp)
679 xp = xx + m_one /
sqrt(m_two)
680 stepf = m_one - m_half *
sqrt(m_e) *
exp(-xp * xp)
692 type(
smear_t),
intent(in) :: this
693 real(real64),
intent(in) :: xx
695 real(real64),
parameter :: maxarg = 200.0_real64
696 real(real64) :: xp, arg, hd, hp, hpm1, aa
705 select case (this%method)
709 if (abs(xx) <= 36.0_real64)
then
710 xp = m_one / (m_one +
exp(-xx))
711 entropyf = xp *
log(xp) + (m_one - xp) *
log(m_one - xp)
715 xp = xx - m_one /
sqrt(m_two)
716 arg = min(maxarg, xp**2)
718 entropyf = m_one /
sqrt(m_two * m_pi) * xp *
exp(-arg)
721 arg = min(maxarg, xx**2)
722 entropyf = -m_half *
exp(-arg) /
sqrt(m_pi)
724 if (this%MP_n > 0)
then
728 aa = m_one /
sqrt(m_pi)
730 hd = m_two * xx * hp - m_two * ni * hd
733 hp = m_two * xx * hd - m_two * ni * hp
735 aa = -aa / (m_four * ii)
736 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
741 xp = abs(xx) + m_one /
sqrt(m_two)
742 entropyf = -
sqrt(m_e) * (abs(xx) *
exp(-xp * xp) / m_two +
sqrt(m_pi) / m_four * loct_erfc(xp))
752 type(
smear_t),
intent(in) :: this
759 type(
smear_t),
intent(in) :: this
760 type(namespace_t),
intent(in) :: namespace
761 integer,
optional,
intent(in) :: iunit
765 if (this%photodop)
then
766 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy (valence ) = ", &
767 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
768 write(message(2),
'(a,f12.6,1x,a)')
"Fermi energy (conduction) = ", &
769 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
770 call messages_info(2, iunit=iunit, namespace=namespace)
772 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy = ", &
773 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
774 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 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_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)
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
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
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_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)
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