Algorithms
IntegralProblem algorithms
AutoBZCore.IntegralAlgorithm — TypeIntegralAlgorithmAbstract supertype for integration algorithms.
Quadrature
AutoBZCore.QuadratureFunction — TypeQuadratureFunction(; fun=trapz, npt=50)Quadrature rule for the standard interval [-1,1] computed from a function x, w = fun(npt). The nodes and weights should be set so the integral of f on [-1,1] is sum(w .* f.(x)). The default quadrature rule is trapz, although other packages provide rules, e.g.
using FastGaussQuadrature
alg = QuadratureFunction(fun=gausslegendre, npt=100)AutoBZCore.QuadGKJL — TypeQuadGKJL(; order = 7, norm = norm)Duplicate of the QuadGKJL provided by Integrals.jl.
AutoBZCore.AuxQuadGKJL — TypeAuxQuadGKJL(; order = 7, norm = norm)Generalization of the QuadGKJL provided by Integrals.jl that allows for AuxValued integrands for auxiliary integration and multi-threaded evaluation with the batch argument to IntegralProblem
Cubature
AutoBZCore.HCubatureJL — TypeHCubatureJL(; norm=norm, initdiv=1)Multi-dimensional h-adaptive cubature from HCubature.jl.
AutoBZCore.MonkhorstPack — TypeMonkhorstPack(; npt=50, syms=nothing)Periodic trapezoidal rule with a fixed number of k-points per dimension, npt, using the PTR rule from AutoSymPTR.jl. The caller should check that the integral is converged w.r.t. npt.
AutoBZCore.AutoSymPTRJL — TypeAutoSymPTRJL(; norm=norm, a=1.0, nmin=50, nmax=1000, n₀=6, Δn=log(10), keepmost=2)Periodic trapezoidal rule with automatic convergence to tolerances passed to the solver with respect to norm using the routine autosymptr from AutoSymPTR.jl. This algorithm is the most efficient for smooth integrands.
Meta-algorithms
AutoBZCore.NestedQuad — TypeNestedQuad(alg::IntegralAlgorithm)
NestedQuad(algs::IntegralAlgorithm...)Nested integration by repeating one quadrature algorithm or composing a list of algorithms. The domain of integration must be an AbstractIteratedLimits from the IteratedIntegration.jl package. Analogous to nested_quad from IteratedIntegration.jl. The integrand should expect SVector inputs. Do not use this for very high-dimensional integrals, since the compilation time scales very poorly with respect to dimensionality.
AutoBZCore.EvalCounter — TypeEvalCounter(alg::IntegralAlgorithm)This algorithm wrapper will count the number of function evaluations used by integration algorithm alg during a solve. This information is found in the sol.stats.numevals field of an IntegralSolution.
AutoBZCore.EvalLogger — TypeEvalLogger(alg::IntegralAlgorithm)This algorithm wrapper will record the quadrature points used by integration algorithm alg during a solve. This information is found in the sol.stats.evallog field of an IntegralSolution.
AutoBZProblem algorithms
Although different algorithms may use different representations of the BZ, the BZ loaded from load_bz can be called with any of the algorithms below, which are aliases for algorithms above
AutoBZCore.AutoBZAlgorithm — TypeAutoBZAlgorithmAbstract supertype for Brillouin zone integration algorithms. All integration problems on the BZ get rescaled to fractional coordinates so that the Brillouin zone becomes [0,1]^d, and integrands should have this periodicity. If the integrand depends on the Brillouin zone basis, then it may have to be transformed to the Cartesian coordinates as a post-processing step. These algorithms also use the symmetries of the Brillouin zone and the integrand.
AutoBZCore.IAI — TypeIAI(alg::IntegralAlgorithm=AuxQuadGKJL())
IAI(algs::IntegralAlgorithm...)Iterated-adaptive integration using nested_quad from IteratedIntegration.jl. This algorithm is the most efficient for localized integrands.
AutoBZCore.TAI — TypeTAI(; norm=norm, initdivs=1)Tree-adaptive integration using hcubature from HCubature.jl. This routine is limited to integration over hypercube domains and may not use all symmetries.
AutoBZCore.PTR — TypePTR(; npt=50)Periodic trapezoidal rule with a fixed number of k-points per dimension, npt, using the routine ptr from AutoSymPTR.jl. The caller should check that the integral is converged w.r.t. npt.
AutoBZCore.AutoPTR — TypeAutoPTR(; norm=norm, a=1.0, nmin=50, nmax=1000, n₀=6, Δn=log(10), keepmost=2)Periodic trapezoidal rule with automatic convergence to tolerances passed to the solver with respect to norm using the routine autosymptr from AutoSymPTR.jl. This algorithm is the most efficient for smooth integrands.
DOSProblem algorithms
Currently the available algorithms are an initial release and we would like to include the following reference algorithms that are also common in the literature in a future release:
- Adaptive Gaussian broadening
AutoBZCore.DOSAlgorithm — TypeDOSAlgorithmAbstract supertype for algorithms for computing density of states
AutoBZCore.GGR — TypeGGR(; npt=50)Generalized Gilat-Raubenheimer method as in "Generalized Gilat–Raubenheimer method for density-of-states calculation in photonic crystals". This method requires the Hamiltonian and its derivatives, and performs a linear extrapolation at each k-point in an equispace grid. The algorithm is expected to show second-order convergence and suffer reduced error at band crossings compared to interpolatory methods.
Arguments
npt: the number of k-points per dimension
AutoBZCore.BCD — TypeBCD(npt, α, ΔE)Brillouin Contour Deformation method was developed in "Efficient extraction of resonant states in systems with defects". This method requires the Hamiltonian and its first and second order derivatives. It performs a deformation of the Brillouin zone into the complex planes based on the first derivative of the Hamiltonian at singularity points. The deformation is controlled by two parameters α and ΔE which scale the deformation around singularity respectively in amplitude and width. This method is expected to show exponential convergence and is not highly sensitive to the choice of parameters. Therefore it can be used with default parameters in most cases.
Arguments
npt: the number of k-points per dimensionα: a scaling parameter for the amplitude of the deformationΔE: a parameter which impacts the width of the deformation at singularities
AutoBZCore.ImplicitIntegrationJL — TypeImplicitIntegrationJL(; kws...)This algorithm is implemented in an extension. Try it with using ImplicitIntegration.
AutoBZCore.LT — TypeLT(npt)Linear Tetrahedron method "High-precision sampling for Brillouin-zone integration in metals". This method requires Hamiltonian's eigenvalues. It performs a linear interpolation of the eigenvalues on a tetrahedric decomposition of the Brillouin zone. Therefore only the eigenvalues are needed at each k-point to perform interpolation. This method is expected to show quadratic convergence.
Arguments
npt: the number of k-points per dimension