Integrands

AutoBZ.jl defines integrands to compute various physical observables, including the density of states, electronic density, and optical conductivity. To define new observables, visit AutoBZCore.jl for general-purpose interfaces to define integrals.

Constructors

AutoBZ.GlocSolverFunction
GlocSolver(Σ, h, bz, bzalg, [linalg=JLInv()]; ω, μ=0, kws...)

Green's function integrands accepting a self energy Σ that can either be a matrix or a function of ω (see the self energy section of the documentation for examples) Additional keywords are passed directly to the solver. Use AutoBZ.update_gloc!(solver; ω, μ=0) to change the parameters. The linalg argument sets the linear system solver

source
AutoBZ.TrGlocSolverFunction
TrGlocSolver(Σ, h, bz, bzalg, [trinvalg=JLTrInv()]; ω, μ=0, kws...)

Green's function integrands accepting a self energy Σ that can either be a matrix or a function of ω (see the self energy section of the documentation for examples) Additional keywords are passed directly to the solver. Use AutoBZ.update_trgloc!(solver; ω, μ=0) to change the parameters.

source
AutoBZ.DOSSolverFunction
DOSSolver(Σ, h, bz, bzalg, [trinvalg=JLTrInv()]; ω, μ=0, kws...)

Green's function integrands accepting a self energy Σ that can either be a matrix or a function of ω (see the self energy section of the documentation for examples) Additional keywords are passed directly to the solver. Use AutoBZ.update_dos!(solver; ω, μ=0) to change the parameters.

source
AutoBZ.TransportFunctionSolverFunction
TransportFunctionSolver(hv::AbstractVelocityInterp, bz, bzalg; β, μ=0, kernel=:fermi, kws...)

A function whose integral over the BZ gives the transport function, proportional to the Drude weight,

\[D_{\alpha\beta} = \sum_{nm} \int_{\text{BZ}} dk f'(\epsilon_{nk}-\mu) \nu_{n\alpha}(k) \nu_{m\beta}(k)\]

where $f(\omega) = (e^{\beta\omega}+1)^{-1}$ is the Fermi distribution. Additional keywords are passed directly to the solver. Use AutoBZ.update_tf!(solver; β, μ=0) to update the parameters.

If the keyword kernel is set to :lorentzian then the following is computed

\[D_{\alpha\beta} = \sum_{nm} \int_{\text{BZ}} dk \operatorname{Tr}[\nu_{n\alpha}(k) A(k, \mu) \nu_{m\beta}(k)]\]

source
AutoBZ.TransportDistributionSolverFunction
TransportDistributionSolver(Σ, hv::AbstractVelocityInterp, bz, bzalg; ω₁, ω₂, μ=0, kws...)

A function whose integral over the BZ gives the transport distribution

\[\Gamma_{\alpha\beta}(\omega_1, \omega_2) = \int_{\text{BZ}} dk \operatorname{Tr}[\nu_\alpha(k) A(k,\omega_1) \nu_\beta(k) A(k, \omega_2)]\]

Based on TRIQS. Additional keywords are passed directly to the solver. Use AutoBZ.update_td!(solver; ω₁, ω₂, μ=0) to update the parameters.

source
AutoBZ.KineticCoefficientSolverFunction
KineticCoefficientSolver(hv, bz, bzalg, Σ, [fdom,] falg, [linalg=JLInv()]; n, β, Ω, μ=0, scale_inner=inv(abs(det(bz.B))*nsyms(bz)), kws...)
KineticCoefficientSolver(Σ, [fdom,] falg, hv, bz, bzalg, [linalg=JLInv()]; n, β, Ω, μ=0, scale_inner=1, kws...)

A solver for kinetic coefficients. The two orderings of arguments correspond to orders of integration. (The outer integral appears first in the argument list.) Use AutoBZ.update_kc!(solver; β, Ω, μ, n) to change parameters. linalg selects the algorithm to compute the resolvent.

Mathematically, this computes

\[A_{n,\alpha\beta}(\Omega) = \int_{-\infty}^{\infty} d \omega (\beta\omega)^{n} \frac{f(\omega) - f(\omega+\Omega)}{\Omega} \Gamma_{\alpha\beta}(\omega, \omega+\Omega)\]

where $f(\omega) = (e^{\beta\omega}+1)^{-1}$ is the Fermi distriubtion. Based on TRIQS.

source
AutoBZ.OpticalConductivitySolverFunction
OpticalConductivitySolver(hv, bz, bzalg, Σ, [fdom,] falg, [linalg=JLInv()]; β, Ω, μ=0, scale_inner=inv(abs(det(bz.B))*nsyms(bz)), kws...)
OpticalConductivitySolver(Σ, [fdom,] falg, hv, bz, bzalg, [linalg=JLInv()]; β, Ω, μ=0, scale_inner=1, kws...)

A solver for the optical conductivity. For details see KineticCoefficientSolver and note that by default the parameter n=0. Use AutoBZ.update_oc!(solver; β, Ω, μ) to change parameters.

source
AutoBZ.ElectronDensitySolverFunction
ElectronDensitySolver(h, bz, bzalg, Σ, [fdom,] falg, [trinvalg=JLTrInv()]; β, μ=0, scale_inner=inv(abs(det(bz.B))*nsyms(bz)), kws...)
ElectronDensitySolver(Σ, [fdom,] falg, h, bz, bzalg, [trinvalg=JLTrInv()]; β, μ=0, scale_inner=inv(oneunit(μ)), kws...)

A solver for the electron density. The two orderings of arguments correspond to orders of integration. (The outer integral appears first in the argument list.) Use AutoBZ.update_density!(solver; β, μ=0). If fdom is not specified the default is (AutoBZ.lb(Σ), AutoBZ.ub(Σ)). The scale_inner keyword is used to rescale inner integration tolerances and can be estimated to reduce computational effort, although the default will be robust. trinvalg may be specified as an algorithm for the inverse trace calculation.

Mathematically, this computes the electron density:

\[n(\mu) = \int_{-\infty}^{\infty} d \omega f(\omega) \operatorname{DOS}(\omega+\mu)\]

where $f(\omega) = (e^{\beta\omega}+1)^{-1}$ is the Fermi distriubtion. To get the density/number of electrons, multiply the result of this integral by n_sp/det(bz.B)

source
AutoBZ.AuxTransportDistributionSolverFunction
AuxTransportDistributionSolver([auxfun], Σ, hv::AbstractVelocityInterp, bz, bzalg; ω₁, ω₂, μ=0, kws...)

A function whose integral over the BZ gives the transport distribution

\[\Gamma_{\alpha\beta}(\omega_1, \omega_2) = \int_{\text{BZ}} dk \operatorname{Tr}[\nu_\alpha(k) A(k,\omega_1) \nu_\beta(k) A(k, \omega_2)]\]

Based on TRIQS. Additional keywords are passed directly to the solver. Use AutoBZ.update_auxtd!(solver; ω₁, ω₂, μ) to update the parameters.

source
AutoBZ.AuxKineticCoefficientSolverFunction
AuxKineticCoefficientSolver([auxfun], hv, bz, bzalg, Σ, [fdom,] falg, [linalg=JLInv()]; n, β, Ω, μ=0, scale_inner=inv(abs(det(bz.B))*nsyms(bz)), kws...)
AuxKineticCoefficientSolver([auxfun], Σ, [fdom,] falg, hv, bz, bzalg, [linalg=JLInv()]; n, β, Ω, μ=0, scale_inner=1, kws...)

A solver for kinetic coefficients using an auxiliary integrand. The two orderings of arguments correspond to orders of integration. (The outer integral appears first in the argument list.) The default auxfun is the sum of the Green's functions. Use AutoBZ.update_auxkc!(solver; β, Ω, μ, n) to change parameters. If fdom is not specified the default is (AutoBZ.lb(Σ), AutoBZ.ub(Σ)).

source
AutoBZ.AuxOpticalConductivitySolverFunction
AuxOpticalConductivitySolver([auxfun], hv, bz, bzalg, Σ, [fdom,] falg, [linalg=JLInv()]; β, Ω, μ=0, scale_inner=inv(abs(det(bz.B))*nsyms(bz)), kws...)
AuxOpticalConductivitySolver([auxfun], Σ, [fdom,] falg, hv, bz, bzalg, [linalg=JLInv()]; β, Ω, μ=0, scale_inner=1, kws...)

A solver for the optical conductivity. For details see AuxKineticCoefficientSolver and note that by default the parameter n=0. Use AutoBZ.update_auxoc!(solver; β, Ω, μ) to change parameters.

source

Functions

AutoBZ.fermiFunction
fermi(β, ω)
fermi(x)

Evaluates a Fermi distribution with unitless input

\[f(x) = \frac{1}{e^{x}+1}\]

source
AutoBZ.fermi′Function
fermi′(x)

Evaluates a first derivative of the Fermi distribution with unitless input

\[\partial_{x} f(x) = -\frac{1}{2(\cosh(x)+1)}\]

Note that the analytic expression above can be rewritten many ways using hypertrigonometric identities.

source
AutoBZ.fermi_windowFunction
fermi_window(β, ω, Ω)
fermi_window(x, y)

Evaluates a unitless window function with unitless inputs determined by the Fermi distribution $f$ and defined by

\[\chi(x, y) = \frac{f(x) - f(x+y)}{y}\]

In the case y==0 then this simplifies to the derivative of the Fermi distribution.

See also fermi and fermi′.

source
AutoBZ.fermi_window_limitsFunction
fermi_window_limits(Ω, β [; atol=0.0, rtol=1e-20])

Returns limits (a,b) over ω restricted to the interval where the Fermi window is larger than max(atol,rtol*fermi_window(0,β*Ω)). Choosing atol and rtol wisely is important to integrating the entire region of interest, since this is a truncation of an infinite interval, and should be tested for convergence.

source