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
250 qtot, nik, nst, kweights)
251 type(
smear_t),
intent(inout) :: this
253 real(real64),
intent(in) :: eigenvalues(:,:), occupations(:,:)
254 real(real64),
intent(in) :: qtot, kweights(:)
255 integer,
intent(in) :: nik, nst
257 integer,
parameter :: nitmax = 200
258 real(real64),
parameter :: tol = 1.0e-10_real64
259 integer :: ist, ik, iter, maxq, weight, sumq_int
260 real(real64) :: sumq_frac, sum_weight
262 real(real64),
allocatable :: eigenval_list(:)
263 integer,
allocatable :: k_list(:), reorder(:)
264 integer :: fermi_count_up, fermi_count_down
268 maxq = this%el_per_state * nst * this%nspins
269 if (maxq - qtot <= -tol)
then
270 message(1) =
'Not enough states'
271 write(
message(2),
'(6x,a,f12.6,a,i10)')
'(total charge = ', qtot, &
272 ' max charge = ', maxq
278 ist_cycle:
do ist = nst, 1, -1
280 if (occupations(ist, ik) > 1e-5_real64)
then
281 this%e_fermi = eigenvalues(ist, ik)
282 this%ef_occ = occupations(ist, ik) / this%el_per_state
289 else if (this%method == smear_semiconductor)
then
291 safe_allocate(eigenval_list(1:nst * nik))
292 safe_allocate( k_list(1:nst * nik))
293 safe_allocate( reorder(1:nst * nik))
298 eigenval_list(iter) = eigenvalues(ist, ik)
305 call sort(eigenval_list, reorder)
307 sumq_int = int(qtot) * this%nik_factor
308 sumq_frac = qtot * this%nik_factor - sumq_int
310 do iter = 1, nst * nik
311 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor +
m_half)
312 if (.not. weight > 0) cycle
313 this%e_fermi = eigenval_list(iter)
314 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
316 if (sumq_int - weight * this%el_per_state <= 0)
then
324 if (iter - fermi_count_down < 1)
exit
325 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear)
exit
326 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor +
m_half)
328 sumq_int = sumq_int + weight * this%el_per_state
329 sum_weight = sum_weight + weight
331 fermi_count_down = fermi_count_down + 1
334 if (iter + fermi_count_up > nst*nik)
exit
335 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear)
exit
336 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor +
m_half)
338 sum_weight = sum_weight + weight
340 fermi_count_up = fermi_count_up + 1
342 this%fermi_count = fermi_count_up + fermi_count_down - 1
343 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
347 sumq_int = sumq_int - weight * this%el_per_state
351 safe_deallocate_a(eigenval_list)
352 safe_deallocate_a(k_list)
353 safe_deallocate_a(reorder)
356 if (this%photodop)
then
359 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
360 this%e_fermi, this%ef_occ)
364 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
365 this%e_fermi_cond, this%ef_occ_cond)
368 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
377 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
378 type(
smear_t),
intent(inout) :: this
380 real(real64),
intent(in) :: dsmear_in
381 real(real64),
intent(in) :: tol
382 integer,
intent(in) :: nik
383 real(real64),
intent(in) :: eigenvalues(:,:)
384 real(real64),
intent(in) :: kweights(:), q_in
385 integer,
intent(in) :: start_band, end_band
386 real(real64),
intent(out) :: e_fermi, ef_occ
388 integer,
parameter :: nitmax = 200
389 integer :: ist, ik, iter
390 real(real64) :: drange, dsmear, emin, emax, xx, sumq
395 dsmear = max(1e-14_real64, dsmear_in)
396 drange = dsmear *
sqrt(-
log(tol * 0.01_real64))
398 emin = minval(eigenvalues) - drange
399 emax = maxval(eigenvalues) + drange
402 e_fermi =
m_half * (emin + emax)
406 do ist = start_band, end_band
407 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
408 sumq = sumq + kweights(ik) * this%el_per_state * &
413 conv = (abs(sumq - q_in) <= tol)
416 if (sumq <= q_in) emin = e_fermi
417 if (sumq >= q_in) emax = e_fermi
423 message(1) =
'Fermi: did not converge.'
432 type(
smear_t),
intent(in) :: this
433 real(real64),
intent(in) :: eigenvalues(:,:)
434 real(real64),
intent(inout) :: occupations(:,:)
435 integer,
intent(in) :: nik, nst
437 integer :: ik, ist, ifermi
438 real(real64) :: dsmear, xx, dsmear_cond
444 else if (this%method == smear_semiconductor)
then
445 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
450 xx = eigenvalues(ist, ik) - this%e_fermi
451 if (xx < -tol_smear)
then
452 occupations(ist, ik) = this%el_per_state
453 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count)
then
454 occupations(ist, ik) = this%ef_occ * this%el_per_state
457 occupations(ist, ik) =
m_zero
464 dsmear = max(1e-14_real64, this%dsmear)
465 if (this%photodop)
then
466 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
469 if (ist < this%photodop_bandmin)
then
471 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
474 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
482 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
494 real(real64) function smear_calc_entropy(this, eigenvalues, &
495 nik, nst, kweights, occ) result(entropy)
496 type(
smear_t),
intent(inout) :: this
497 real(real64),
intent(in) :: eigenvalues(:,:)
498 real(real64),
intent(in) :: kweights(:)
499 integer,
intent(in) :: nik, nst
500 real(real64),
intent(in) :: occ(:, :)
503 real(real64) :: dsmear, xx, term, ff
505 push_sub(smear_calc_entropy)
509 if (this%method ==
smear_fixed_occ .or. this%method == smear_semiconductor)
then
515 ff = occ(ist, ik) / this%el_per_state
517 term = ff *
log(ff) + (1 - ff) *
log(1 - ff)
521 entropy = entropy - kweights(ik) * this%el_per_state * term
525 dsmear = max(1e-14_real64, this%dsmear)
529 if (eigenvalues(ist, ik) < huge(
m_one))
then
530 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
531 entropy = entropy - kweights(ik) * this%el_per_state * &
538 pop_sub(smear_calc_entropy)
544 type(
smear_t),
intent(in) :: this
545 real(real64),
intent(in) :: xx
547 real(real64),
parameter :: maxarg = 200.0_real64
548 real(real64) :: xp, arg, hd, hp, aa
557 select case (this%method)
559 if (abs(xx) <= m_epsilon)
then
564 if (abs(xx) <= 36.0_real64)
then
565 deltaf = m_one / (m_two +
exp(-xx) +
exp(xx))
569 xp = xx - m_one /
sqrt(m_two)
570 arg = min(maxarg, xp**2)
572 deltaf =
exp(-arg) /
sqrt(m_pi) * (m_two -
sqrt(m_two) * xx)
575 arg = min(maxarg, xx**2)
576 deltaf =
exp(-arg) /
sqrt(m_pi)
578 if (this%MP_n > 0)
then
582 aa = m_one /
sqrt(m_pi)
584 hd = m_two * xx * hp - m_two * ni * hd
586 aa = -aa / (m_four * ii)
587 hp = m_two * xx * hd - m_two * ni * hp
589 deltaf = deltaf + aa * hp
594 xp = abs(xx) + m_one /
sqrt(m_two)
595 deltaf =
sqrt(m_e) * xp *
exp(-xp * xp)
604 type(
smear_t),
intent(in) :: this
605 real(real64),
intent(in) :: xx
607 real(real64),
parameter :: maxarg = 200.0_real64
608 real(real64) :: xp, arg, hd, hp, aa
617 select case (this%method)
619 if (xx > m_zero)
then
621 else if (abs(xx) <= m_epsilon)
then
626 if (xx > maxarg)
then
628 else if (xx > -maxarg)
then
629 stepf = m_one / (m_one +
exp(-xx))
633 xp = xx - m_one /
sqrt(m_two)
634 arg = min(maxarg, xp**2)
636 stepf = m_half * loct_erf(xp) + &
637 m_one /
sqrt(m_two * m_pi) *
exp(-arg) + m_half
640 stepf = m_half * loct_erfc(-xx)
642 if (this%MP_n > 0)
then
644 arg = min(maxarg, xx**2)
647 aa = m_one /
sqrt(m_pi)
649 hd = m_two * xx * hp - m_two * ni * hd
651 aa = -aa / (m_four * ii)
652 stepf = stepf - aa * hd
653 hp = m_two * xx * hd - m_two * ni * hp
659 if (xx <= m_zero)
then
660 xp = xx - m_one /
sqrt(m_two)
661 stepf = m_half *
sqrt(m_e) *
exp(-xp * xp)
663 xp = xx + m_one /
sqrt(m_two)
664 stepf = m_one - m_half *
sqrt(m_e) *
exp(-xp * xp)
676 type(
smear_t),
intent(in) :: this
677 real(real64),
intent(in) :: xx
679 real(real64),
parameter :: maxarg = 200.0_real64
680 real(real64) :: xp, arg, hd, hp, hpm1, aa
689 select case (this%method)
693 if (abs(xx) <= 36.0_real64)
then
694 xp = m_one / (m_one +
exp(-xx))
695 entropyf = xp *
log(xp) + (m_one - xp) *
log(m_one - xp)
699 xp = xx - m_one /
sqrt(m_two)
700 arg = min(maxarg, xp**2)
702 entropyf = m_one /
sqrt(m_two * m_pi) * xp *
exp(-arg)
705 arg = min(maxarg, xx**2)
706 entropyf = -m_half *
exp(-arg) /
sqrt(m_pi)
708 if (this%MP_n > 0)
then
712 aa = m_one /
sqrt(m_pi)
714 hd = m_two * xx * hp - m_two * ni * hd
717 hp = m_two * xx * hd - m_two * ni * hp
719 aa = -aa / (m_four * ii)
720 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
725 xp = abs(xx) + m_one /
sqrt(m_two)
726 entropyf = -
sqrt(m_e) * (abs(xx) *
exp(-xp * xp) / m_two +
sqrt(m_pi) / m_four * loct_erfc(xp))
736 type(
smear_t),
intent(in) :: this
743 type(
smear_t),
intent(in) :: this
744 type(namespace_t),
intent(in) :: namespace
745 integer,
optional,
intent(in) :: iunit
749 if (this%photodop)
then
750 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy (valence ) = ", &
751 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
752 write(message(2),
'(a,f12.6,1x,a)')
"Fermi energy (conduction) = ", &
753 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
754 call messages_info(2, iunit=iunit, namespace=namespace)
756 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy = ", &
757 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
758 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__
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