Integrands
The design of AutoBZCore.jl uses multiple dispatch to provide multiple interfaces for user integrands that allow various optimizations to be compatible with a common interface for solvers.
AutoBZCore.IntegralFunction — TypeIntegralFunction(f, [prototype=nothing, executor=SerialExecutor()])Constructor for an out-of-place integrand of the form f(x, p). Optionally, a prototype can be provided for the output of the function.
AutoBZCore.InplaceIntegralFunction — TypeInplaceIntegralFunction(f!, prototype::AbstractArray)Constructor for an inplace integrand of the form f!(y, x, p). A prototype array is required to store the same type and size as the result, y.
AutoBZCore.InplaceBatchIntegralFunction — TypeInplaceBatchIntegralFunction(f!, prototype; max_batch::Integer=typemax(Int))Constructor for an inplace, batched integrand of the form f!(y, x, p) that accepts an array x containing a batch of evaluation points stored along the last axis of the array. A prototype array is required to store the same type and size as the result, y, however the last axis, which is reserved for batching, should contain at least one element. The max_batch keyword sets a soft limit on the number of points batched simultaneously.
AutoBZCore.CommonSolveIntegralFunction — TypeCommonSolveIntegralFunction(solve!, prob, alg, [prototype, specialize, executor]; kws...)Constructor for an integrand that solves a problem defined with the CommonSolve.jl interface, prob, which is instantiated using init(prob, alg; kws...). The solution = solve!(solver, x, p) function supplied by the caller must do the work of the problem, although it need not be a method of CommonSolve.solve! as the out-of-place semantics of passing the arguments x, p can provide a speedup. The prototype argument can help control how much to specialize on the solution type of the problem. By default, specialize=DefaultSpecialize() uses Julia's default heuristics, which can give up on inference in complicated codes. Additionally, FullSpecialize() can obtain the fastest run times with the longest compile times, NoSpecialize() strikes a good balance of run time, compile time and inference, and FunctionWrapperSpecialize() may have the fastest compile time and very good run times (comparable to FullSpecialize()) but with possible issues regarding world age. The executor keyword specifies how to schedule and run the integrand evaluation, defaulting to SerialExecutor() with an additional option for ThreadedExecutor(::Integer).
AutoBZCore.FourierIntegralFunction — TypeFourierIntegralFunction(f, s, [prototype=nothing, executor=SerialExecutor()]; alias=false)Arguments
f: The integrand, accepting inputsf(x, s(x), p)s::AbstractFourierSeries: The Fourier series to evaluateprototype:alias::Bool: whether todeepcopythe series (false) or use the series as-is (true)
AutoBZCore.CommonSolveFourierIntegralFunction — TypeCommonSolveFourierIntegralFunction(solve!, prob, alg, s, [prototype, specialize, executor]; alias=false, kws...)Constructor for an integrand that solves a problem defined with the CommonSolve.jl interface, prob, which is instantiated using init(prob, alg; kws...). The function sol = solve!(solver, x, s(x), p) should return the value of the solution and s must be a Fourier series. The prototype argument can help control how much to specialize on the solution type of the problem. By default, specialize=DefaultSpecialize() uses Julia's default heuristics, which can give up on inference in complicated codes. Additionally, FullSpecialize() can obtain the fastest run times with the longest compile times, NoSpecialize() strikes a good balance of run time, compile time and inference, and FunctionWrapperSpecialize() may have the fastest compile time and very good run times (comparable to FullSpecialize()) but with possible issues regarding world age. The executor keyword specifies how to schedule and run the integrand evaluation, defaulting to SerialExecutor() with an additional option for ThreadedExecutor(::Integer).