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
66 integer,
public :: MP_n
67 integer :: fermi_count
72 logical,
public:: photodop
73 real(real64) :: nephotodop
74 integer :: photodop_bandmin
78 integer,
parameter,
public :: &
79 SMEAR_SEMICONDUCTOR = 1, &
88 real(real64),
parameter :: TOL_SMEAR = 1e-6_real64
93 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
94 type(smear_t),
intent(out) :: this
95 type(namespace_t),
intent(in) :: namespace
96 integer,
intent(in) :: ispin
97 logical,
intent(in) :: fixed_occ
98 logical,
intent(in) :: integral_occs
99 type(kpoints_t),
intent(in) :: kpoints
103 this%integral_occs = integral_occs
135 call parse_variable(namespace,
'SmearingFunction', smear_semiconductor, this%method)
149 this%dsmear = 1e-14_real64
150 if (this%method /= smear_semiconductor .and. this%method /=
smear_fixed_occ)
then
163 this%dsmear_cond = this%dsmear
181 this%photodop=.false.
192 call parse_variable(namespace,
'PhotoDopingBand', 1, this%photodop_bandmin)
198 this%el_per_state = 2
200 this%el_per_state = 1
209 if (this%method == smear_semiconductor)
then
212 if (this%nik_factor == 0)
then
213 message(1) =
"k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
236 type(
smear_t),
intent(out) :: to
237 type(
smear_t),
intent(in) :: from
241 to%method = from%method
242 to%dsmear = from%dsmear
243 to%e_fermi = from%e_fermi
244 to%el_per_state = from%el_per_state
245 to%ef_occ = from%ef_occ
247 to%fermi_count = from%fermi_count
248 to%nik_factor = from%nik_factor
256 qtot, nik, nst, kweights)
257 type(
smear_t),
intent(inout) :: this
259 real(real64),
intent(in) :: eigenvalues(:,:), occupations(:,:)
260 real(real64),
intent(in) :: qtot, kweights(:)
261 integer,
intent(in) :: nik, nst
263 real(real64),
parameter :: tol = 1.0e-10_real64
264 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
265 real(real64) :: sumq_frac
267 real(real64),
allocatable :: eigenval_list(:)
268 integer,
allocatable :: k_list(:), reorder(:)
269 integer :: fermi_count_up, fermi_count_down
273 maxq = this%el_per_state * nst * this%nspins
274 if (maxq - qtot <= -tol)
then
275 message(1) =
'Not enough states'
276 write(
message(2),
'(6x,a,f12.6,a,i10)')
'(total charge = ', qtot, &
277 ' max charge = ', maxq
283 ist_cycle:
do ist = nst, 1, -1
284 if (any(occupations(ist, :) >
m_min_occ))
then
285 this%e_fermi = eigenvalues(ist, 1)
286 this%ef_occ = occupations(ist, 1) / this%el_per_state
288 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) >
m_min_occ)
then
289 this%e_fermi = eigenvalues(ist, ik)
290 this%ef_occ = occupations(ist, ik) / this%el_per_state
297 else if (this%method == smear_semiconductor)
then
299 safe_allocate(eigenval_list(1:nst * nik))
300 safe_allocate( k_list(1:nst * nik))
301 safe_allocate( reorder(1:nst * nik))
306 eigenval_list(iter) = eigenvalues(ist, ik)
313 call sort(eigenval_list, reorder)
315 sumq_int =
floor(qtot) * this%nik_factor
316 sumq_frac = qtot * this%nik_factor - sumq_int
317 if ( sumq_frac + tol > this%nik_factor )
then
318 sumq_int = sumq_int + this%nik_factor
319 sumq_frac = sumq_frac - this%nik_factor
321 if (sumq_frac < tol) sumq_frac =
m_zero
323 do iter = 1, nst * nik
324 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor +
m_half)
325 if (.not. weight > 0) cycle
326 this%e_fermi = eigenval_list(iter)
327 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
329 if (sumq_int - weight * this%el_per_state <= 0)
then
337 if (iter - fermi_count_down < 1)
exit
338 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear)
exit
339 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor +
m_half)
341 sumq_int = sumq_int + weight * this%el_per_state
342 sum_weight = sum_weight + weight
344 fermi_count_down = fermi_count_down + 1
347 if (iter + fermi_count_up > nst*nik)
exit
348 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear)
exit
349 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor +
m_half)
351 sum_weight = sum_weight + weight
353 fermi_count_up = fermi_count_up + 1
355 this%fermi_count = fermi_count_up + fermi_count_down - 1
356 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
360 sumq_int = sumq_int - weight * this%el_per_state
364 safe_deallocate_a(eigenval_list)
365 safe_deallocate_a(k_list)
366 safe_deallocate_a(reorder)
369 if (this%photodop)
then
372 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
373 this%e_fermi, this%ef_occ)
377 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
378 this%e_fermi_cond, this%ef_occ_cond)
381 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
390 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
391 type(
smear_t),
intent(inout) :: this
393 real(real64),
intent(in) :: dsmear_in
394 real(real64),
intent(in) :: tol
395 integer,
intent(in) :: nik
396 real(real64),
intent(in) :: eigenvalues(:,:)
397 real(real64),
intent(in) :: kweights(:), q_in
398 integer,
intent(in) :: start_band, end_band
399 real(real64),
intent(out) :: e_fermi, ef_occ
401 integer,
parameter :: nitmax = 200
402 integer :: ist, ik, iter
403 real(real64) :: drange, dsmear, emin, emax, xx, sumq
408 dsmear = max(1e-14_real64, dsmear_in)
409 drange = dsmear *
sqrt(-
log(tol * 0.01_real64))
411 emin = minval(eigenvalues) - drange
412 emax = maxval(eigenvalues) + drange
415 e_fermi =
m_half * (emin + emax)
419 do ist = start_band, end_band
420 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
421 sumq = sumq + kweights(ik) * this%el_per_state * &
426 conv = (abs(sumq - q_in) <= tol)
429 if (sumq <= q_in) emin = e_fermi
430 if (sumq >= q_in) emax = e_fermi
436 message(1) =
'Fermi: did not converge.'
445 type(
smear_t),
intent(in) :: this
446 real(real64),
intent(in) :: eigenvalues(:,:)
447 real(real64),
intent(inout) :: occupations(:,:)
448 integer,
intent(in) :: nik, nst
450 integer :: ik, ist, ifermi
451 real(real64) :: dsmear, xx, dsmear_cond
457 else if (this%method == smear_semiconductor)
then
458 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
463 xx = eigenvalues(ist, ik) - this%e_fermi
464 if (xx < -tol_smear)
then
465 occupations(ist, ik) = this%el_per_state
466 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count)
then
467 occupations(ist, ik) = this%ef_occ * this%el_per_state
470 occupations(ist, ik) =
m_zero
477 dsmear = max(1e-14_real64, this%dsmear)
478 if (this%photodop)
then
479 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
482 if (ist < this%photodop_bandmin)
then
484 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
487 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
495 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
507 real(real64) function smear_calc_entropy(this, eigenvalues, &
508 nik, nst, kweights, occ) result(entropy)
509 type(
smear_t),
intent(inout) :: this
510 real(real64),
intent(in) :: eigenvalues(:,:)
511 real(real64),
intent(in) :: kweights(:)
512 integer,
intent(in) :: nik, nst
513 real(real64),
intent(in) :: occ(:, :)
516 real(real64) :: dsmear, xx, term, ff
518 push_sub(smear_calc_entropy)
522 if (this%method ==
smear_fixed_occ .or. this%method == smear_semiconductor)
then
528 ff = occ(ist, ik) / this%el_per_state
530 term = ff *
log(ff) + (1 - ff) *
log(1 - ff)
534 entropy = entropy - kweights(ik) * this%el_per_state * term
538 dsmear = max(1e-14_real64, this%dsmear)
542 if (eigenvalues(ist, ik) < huge(
m_one))
then
543 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
544 entropy = entropy - kweights(ik) * this%el_per_state * &
551 pop_sub(smear_calc_entropy)
557 type(
smear_t),
intent(in) :: this
558 real(real64),
intent(in) :: xx
560 real(real64),
parameter :: maxarg = 200.0_real64
561 real(real64) :: xp, arg, hd, hp, aa
570 select case (this%method)
572 if (abs(xx) <= m_epsilon)
then
577 if (abs(xx) <= 36.0_real64)
then
578 deltaf = m_one / (m_two +
exp(-xx) +
exp(xx))
582 xp = xx - m_one /
sqrt(m_two)
583 arg = min(maxarg, xp**2)
585 deltaf =
exp(-arg) /
sqrt(m_pi) * (m_two -
sqrt(m_two) * xx)
588 arg = min(maxarg, xx**2)
589 deltaf =
exp(-arg) /
sqrt(m_pi)
591 if (this%MP_n > 0)
then
595 aa = m_one /
sqrt(m_pi)
597 hd = m_two * xx * hp - m_two * ni * hd
599 aa = -aa / (m_four * ii)
600 hp = m_two * xx * hd - m_two * ni * hp
602 deltaf = deltaf + aa * hp
607 xp = abs(xx) + m_one /
sqrt(m_two)
608 deltaf =
sqrt(m_e) * xp *
exp(-xp * xp)
611 deltaf =
exp(-xx * xx) /
sqrt(m_pi)
614 deltaf = m_one / (m_pi * (xx * xx + m_one))
623 type(
smear_t),
intent(in) :: this
624 real(real64),
intent(in) :: xx
626 real(real64),
parameter :: maxarg = 200.0_real64
627 real(real64) :: xp, arg, hd, hp, aa
636 select case (this%method)
638 if (xx > m_zero)
then
640 else if (abs(xx) <= m_epsilon)
then
645 if (xx > maxarg)
then
647 else if (xx > -maxarg)
then
648 stepf = m_one / (m_one +
exp(-xx))
652 xp = xx - m_one /
sqrt(m_two)
653 arg = min(maxarg, xp**2)
655 stepf = m_half * erf(xp) + &
656 m_one /
sqrt(m_two * m_pi) *
exp(-arg) + m_half
659 stepf = m_half * erfc(-xx)
661 if (this%MP_n > 0)
then
663 arg = min(maxarg, xx**2)
666 aa = m_one /
sqrt(m_pi)
668 hd = m_two * xx * hp - m_two * ni * hd
670 aa = -aa / (m_four * ii)
671 stepf = stepf - aa * hd
672 hp = m_two * xx * hd - m_two * ni * hp
678 if (xx <= m_zero)
then
679 xp = xx - m_one /
sqrt(m_two)
680 stepf = m_half *
sqrt(m_e) *
exp(-xp * xp)
682 xp = xx + m_one /
sqrt(m_two)
683 stepf = m_one - m_half *
sqrt(m_e) *
exp(-xp * xp)
687 stepf = m_half * erfc(-xx)
690 stepf =
atan(xx) / m_pi + m_half
701 type(
smear_t),
intent(in) :: this
702 real(real64),
intent(in) :: xx
704 real(real64),
parameter :: maxarg = 200.0_real64
705 real(real64) :: xp, arg, hd, hp, hpm1, aa
714 select case (this%method)
718 if (abs(xx) <= 36.0_real64)
then
719 xp = m_one / (m_one +
exp(-xx))
720 entropyf = xp *
log(xp) + (m_one - xp) *
log(m_one - xp)
724 xp = xx - m_one /
sqrt(m_two)
725 arg = min(maxarg, xp**2)
727 entropyf = m_one /
sqrt(m_two * m_pi) * xp *
exp(-arg)
730 arg = min(maxarg, xx**2)
731 entropyf = -m_half *
exp(-arg) /
sqrt(m_pi)
733 if (this%MP_n > 0)
then
737 aa = m_one /
sqrt(m_pi)
739 hd = m_two * xx * hp - m_two * ni * hd
742 hp = m_two * xx * hd - m_two * ni * hp
744 aa = -aa / (m_four * ii)
745 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
750 xp = abs(xx) + m_one /
sqrt(m_two)
751 entropyf = -
sqrt(m_e) * (abs(xx) *
exp(-xp * xp) / m_two +
sqrt(m_pi) / m_four * erfc(xp))
754 entropyf = -
exp(-xx * xx) / m_two /
sqrt(m_pi)
767 type(
smear_t),
intent(in) :: this
774 type(
smear_t),
intent(in) :: this
775 type(namespace_t),
intent(in) :: namespace
776 integer,
optional,
intent(in) :: iunit
780 if (this%photodop)
then
781 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy (valence ) = ", &
782 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
783 write(message(2),
'(a,f12.6,1x,a)')
"Fermi energy (conduction) = ", &
784 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
785 call messages_info(2, iunit=iunit, namespace=namespace)
787 write(message(1),
'(a,f12.6,1x,a)')
"Fermi energy = ", &
788 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
789 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_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
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