![]() |
(git:1155b05)
|
Unified smearing module supporting four methods: smear_fermi_dirac — Fermi-Dirac distribution smear_gaussian — Gaussian broadening smear_mp — Methfessel-Paxton first order smear_mv — Marzari-Vanderbilt (cold smearing) More...
Functions/Subroutines | |
| subroutine, public | smearocc (f, n, kts, e, mu, sigma, maxocc, method, estate, festate) |
| Returns occupations and smearing correction for a given set of energies and chemical potential, using one of four smearing methods. | |
| subroutine, public | smearfixed (f, mu, kts, e, n, sigma, maxocc, method, estate, festate) |
| Bisection search for the chemical potential mu such that the total electron count equals N, for a given smearing method (Gamma point). Brackets mu by expanding outward from [min(e), max(e)] in steps of sigma, then bisects to machine precision. Could fail if mu lies far outside the eigenvalue range (to be fixed). | |
| subroutine, public | smearkp (f, mu, kts, e, nel, wk, sigma, maxocc, method) |
| Bisection search for mu given a target electron count (k-point case, single spin channel or spin-degenerate). Initial bracket width is max(10*sigma, 0.5) for Gaussian/MP/MV, or sigma*ln[(1-eps)/eps] for Fermi-Dirac, reflecting the different tail decay rates. | |
| subroutine, public | smearkp2 (f, mu, kts, e, nel, wk, sigma, method) |
| Bisection search for mu (k-point, spin-polarised with a shared chemical potential across both spin channels). Asserts that the third dimension of f and e is exactly 2. | |
| subroutine, public | smearfixedderiv (dfde, f, mu, kts, e, n, sigma, maxocc, method, estate, festate) |
| Analytical Jacobian df_i/de_j for any smearing method under the electron-number constraint sum(f) = N. | |
| subroutine, public | smearfixedderivmv (result, f, mu, kts, e, n_el, sigma, maxocc, method, v, estate, festate) |
| Apply TRANSPOSE(df/de) to a vector WITHOUT forming the full N x N Jacobian. O(N) time and O(N) memory for all four methods. | |
Unified smearing module supporting four methods: smear_fermi_dirac — Fermi-Dirac distribution smear_gaussian — Gaussian broadening smear_mp — Methfessel-Paxton first order smear_mv — Marzari-Vanderbilt (cold smearing)
All methods share the bisection framework, LFOMO/HOMO logic, and the analytical rank-1 Jacobian. Only the per-state math (f, kTS, g_i) differs, selected via a method integer from input_constants.
| subroutine, public smearing_utils::smearocc | ( | real(kind=dp), dimension(:), intent(out) | f, |
| real(kind=dp), intent(out) | n, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:), intent(in) | e, | ||
| real(kind=dp), intent(in) | mu, | ||
| real(kind=dp), intent(in) | sigma, | ||
| real(kind=dp), intent(in) | maxocc, | ||
| integer, intent(in) | method, | ||
| integer, intent(in), optional | estate, | ||
| real(kind=dp), intent(in), optional | festate | ||
| ) |
Returns occupations and smearing correction for a given set of energies and chemical potential, using one of four smearing methods.
Fermi-Dirac: f_i = occ / [1 + exp((e_i - mu)/sigma)] Gaussian: f_i = (occ/2) * erfc[(e_i - mu)/sigma] MP-1: f_i = (occ/2) * erfc(x) - occ*x/(2*sqrt(pi)) * exp(-x^2) MV: f_i = (occ/2) * erfc(u) + occ/(sqrt(2*pi)) * exp(-u^2), u = x + 1/sqrt(2)
kTS is the smearing correction to the free energy (physically -TS for Fermi-Dirac; a variational correction term for the other methods). It enters the total energy and the Gillan extrapolation E(0) = E - kTS/2.
| f | occupations (output) |
| N | total number of electrons (output) |
| kTS | smearing correction to the free energy (output) |
| e | eigenvalues (input) |
| mu | chemical potential (input) |
| sigma | smearing width: kT for Fermi-Dirac, sigma for others (input) |
| maxocc | maximum occupation of an orbital (input) |
| method | smearing method selector from input_constants (input) |
| estate | excited state index for core-level spectroscopy (optional) |
| festate | occupation of the excited state (optional) |
Definition at line 74 of file smearing_utils.F.
| subroutine, public smearing_utils::smearfixed | ( | real(kind=dp), dimension(:), intent(out) | f, |
| real(kind=dp), intent(out) | mu, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:), intent(in) | e, | ||
| real(kind=dp), intent(in) | n, | ||
| real(kind=dp), intent(in) | sigma, | ||
| real(kind=dp), intent(in) | maxocc, | ||
| integer, intent(in) | method, | ||
| integer, intent(in), optional | estate, | ||
| real(kind=dp), intent(in), optional | festate | ||
| ) |
Bisection search for the chemical potential mu such that the total electron count equals N, for a given smearing method (Gamma point). Brackets mu by expanding outward from [min(e), max(e)] in steps of sigma, then bisects to machine precision. Could fail if mu lies far outside the eigenvalue range (to be fixed).
| f | occupations (output) |
| mu | chemical potential found by bisection (output) |
| kTS | smearing correction (output) |
| e | eigenvalues (input) |
| N | target number of electrons (input) |
| sigma | smearing width (input) |
| maxocc | maximum occupation (input) |
| method | smearing method selector (input) |
| estate | excited state index for core-level spectroscopy (optional) |
| festate | occupation of the excited state (optional) |
Definition at line 276 of file smearing_utils.F.
| subroutine, public smearing_utils::smearkp | ( | real(kind=dp), dimension(:, :), intent(out) | f, |
| real(kind=dp), intent(out) | mu, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:, :), intent(in) | e, | ||
| real(kind=dp), intent(in) | nel, | ||
| real(kind=dp), dimension(:), intent(in) | wk, | ||
| real(kind=dp), intent(in) | sigma, | ||
| real(kind=dp), intent(in) | maxocc, | ||
| integer, intent(in) | method | ||
| ) |
Bisection search for mu given a target electron count (k-point case, single spin channel or spin-degenerate). Initial bracket width is max(10*sigma, 0.5) for Gaussian/MP/MV, or sigma*ln[(1-eps)/eps] for Fermi-Dirac, reflecting the different tail decay rates.
| f | occupations (nmo x nkp, output) |
| mu | chemical potential (output) |
| kTS | smearing correction (output) |
| e | eigenvalues (nmo x nkp, input) |
| nel | target number of electrons (input) |
| wk | k-point weights (input) |
| sigma | smearing width (input) |
| maxocc | maximum occupation (input) |
| method | smearing method selector (input) |
Definition at line 363 of file smearing_utils.F.
| subroutine, public smearing_utils::smearkp2 | ( | real(kind=dp), dimension(:, :, :), intent(out) | f, |
| real(kind=dp), intent(out) | mu, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:, :, :), intent(in) | e, | ||
| real(kind=dp), intent(in) | nel, | ||
| real(kind=dp), dimension(:), intent(in) | wk, | ||
| real(kind=dp), intent(in) | sigma, | ||
| integer, intent(in) | method | ||
| ) |
Bisection search for mu (k-point, spin-polarised with a shared chemical potential across both spin channels). Asserts that the third dimension of f and e is exactly 2.
| f | occupations (nmo x nkp x 2, output) |
| mu | chemical potential (output) |
| kTS | smearing correction (output) |
| e | eigenvalues (nmo x nkp x 2, input) |
| nel | target total number of electrons (input) |
| wk | k-point weights (input) |
| sigma | smearing width (input) |
| method | smearing method selector (input) |
Definition at line 427 of file smearing_utils.F.
| subroutine, public smearing_utils::smearfixedderiv | ( | real(kind=dp), dimension(:, :), intent(out) | dfde, |
| real(kind=dp), dimension(:), intent(out) | f, | ||
| real(kind=dp), intent(out) | mu, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:), intent(in) | e, | ||
| real(kind=dp), intent(in) | n, | ||
| real(kind=dp), intent(in) | sigma, | ||
| real(kind=dp), intent(in) | maxocc, | ||
| integer, intent(in) | method, | ||
| integer, intent(in), optional | estate, | ||
| real(kind=dp), intent(in), optional | festate | ||
| ) |
Analytical Jacobian df_i/de_j for any smearing method under the electron-number constraint sum(f) = N.
Differentiating f_i(e, mu(e)) where mu is implicitly defined by the constraint yields:
df_i/de_j = -delta_{ij} * g_i + g_i * g_j / G
where g_i = -df_i/de_i (mu fixed) and G = sum(g_i). This is a diagonal matrix plus a symmetric rank-1 update. Building it costs O(N) for g, plus O(N^2) for the outer product.
Replaces the original numerical finite-difference FermiFixedDeriv which required 2N bisection solves. Exact to machine precision for all four methods.
| dfde | Jacobian matrix dfde(i,j) = df_i/de_j (Nstate x Nstate, output) |
| f | occupations (output) |
| mu | chemical potential (output) |
| kTS | smearing correction (output) |
| e | eigenvalues (input) |
| N | target number of electrons (input) |
| sigma | smearing width (input) |
| maxocc | maximum occupation (input) |
| method | smearing method selector (input) |
| estate | excited state index (optional) |
| festate | occupation of the excited state (optional) |
Definition at line 592 of file smearing_utils.F.
| subroutine, public smearing_utils::smearfixedderivmv | ( | real(kind=dp), dimension(:), intent(out) | result, |
| real(kind=dp), dimension(:), intent(out) | f, | ||
| real(kind=dp), intent(out) | mu, | ||
| real(kind=dp), intent(out) | kts, | ||
| real(kind=dp), dimension(:), intent(in) | e, | ||
| real(kind=dp), intent(in) | n_el, | ||
| real(kind=dp), intent(in) | sigma, | ||
| real(kind=dp), intent(in) | maxocc, | ||
| integer, intent(in) | method, | ||
| real(kind=dp), dimension(:), intent(in) | v, | ||
| integer, intent(in), optional | estate, | ||
| real(kind=dp), intent(in), optional | festate | ||
| ) |
Apply TRANSPOSE(df/de) to a vector WITHOUT forming the full N x N Jacobian. O(N) time and O(N) memory for all four methods.
Exploiting the rank-1 structure of the constrained Jacobian:
[J^T v]_j = g_j * (g . v / G - v_j)
This replaces the pattern used in qs_mo_occupation: ALLOCATE(dfde(nmo,nmo)) CALL SmearFixedDeriv(dfde, ...) RESULT = MATMUL(TRANSPOSE(dfde), v) DEALLOCATE(dfde) turning O(N^2) storage + O(N^2) MATMUL into O(N) throughout.
Currently the sole caller (qs_ot_scf do_ener) is dead code, but this routine is ready for when it is enabled.
| RESULT | output vector = TRANSPOSE(df/de) * v (Nstate, output) |
| f | occupations (output) |
| mu | chemical potential (output) |
| kTS | smearing correction (output) |
| e | eigenvalues (input) |
| N_el | target number of electrons (input) |
| sigma | smearing width (input) |
| maxocc | maximum occupation (input) |
| method | smearing method selector (input) |
| v | input vector to multiply (Nstate, input) |
| estate | excited state index (optional) |
| festate | occupation of the excited state (optional) |
Definition at line 665 of file smearing_utils.F.