API Documentation
This is the stable public API of BAT. Forward/backward compatibility follows Julia's semantic versioning rules.
Types
BAT.AbstractMCMCWeightingSchemeBAT.AbstractMedianEstimatorBAT.AbstractModeEstimatorBAT.AbstractPosteriorMeasureBAT.AbstractSamplingAlgorithmBAT.AbstractTransformTargetBAT.AdaptiveAffineTuningBAT.AssumeConvergenceBAT.AutocorLenAlgorithmBAT.BATContextBAT.BATHDF5IOBAT.BATIOAlgorithmBAT.BinningAlgorithmBAT.BrooksGelmanConvergenceBAT.ConvergenceTestBAT.CuhreIntegrationBAT.DensitySampleBAT.DensitySampleVectorBAT.DivonneIntegrationBAT.DoNotTransformBAT.EffSampleSizeAlgorithmBAT.EffSampleSizeFromACBAT.EvaluatedMeasureBAT.ExplicitInitBAT.FixedMGVIScheduleBAT.FixedNBinsBAT.FreedmanDiaconisBinningBAT.GelmanRubinConvergenceBAT.GeyerAutocorLenBAT.HamiltonianMCBAT.IIDSamplingBAT.IdentityTransformAlgorithmBAT.InitFromIIDBAT.InitFromSamplesBAT.InitFromTargetBAT.InitvalAlgorithmBAT.IntegrationAlgorithmBAT.KishESSBAT.MCMCAlgorithmBAT.MCMCBurninAlgorithmBAT.MCMCChainPoolInitBAT.MCMCInitAlgorithmBAT.MCMCMultiCycleBurninBAT.MCMCProposalTuningBAT.MCMCRetryInitBAT.MCMCTransformTuningBAT.MGVISamplingBAT.MGVIScheduleBAT.MaxDensitySearchBAT.ModeAsDefinedBAT.NoMCMCProposalTuningBAT.NoMCMCTransformTuningBAT.OptimAlgBAT.OptimizationAlgBAT.OrderedResamplingBAT.PosteriorMeasureBAT.PriorSubstitutionBAT.PriorToNormalBAT.PriorToUniformBAT.RAMTuningBAT.RandResamplingBAT.RandomWalkBAT.RepetitionWeightingBAT.RiceBinningBAT.SampleMedianEstimatorBAT.ScottBinningBAT.SokalAutocorLenBAT.SquareRootBinningBAT.SturgesBinningBAT.SuaveIntegrationBAT.ToRealVectorBAT.TransformAlgorithmBAT.TransformedMCMCBAT.VEGASIntegration
Functions and macros
BAT.bat_convergenceBAT.bat_defaultBAT.bat_eff_sample_sizeBAT.bat_findmedianBAT.bat_findmodeBAT.bat_initvalBAT.bat_integrateBAT.bat_readBAT.bat_sampleBAT.bat_transformBAT.bat_writeBAT.distbindBAT.distprodBAT.get_batcontextBAT.lbqintegralBAT.log_batdebugBAT.set_batcontextBAT.unevaluated
Documentation
BAT.bat_convergence — Functionbat_convergence(
algoutput,
[algorithm::ConvergenceTest],
[context::BATContext]
)Check if an algorithm has converged, based on it's output algoutput
Returns a NamedTuple of the shape
(result, ...)result indicates whether algoutput and must either be a Bool or support convert(Bool, result). It should typically contains measures of algorithm convergence, like a convergence value and it's threshold, etc.
Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_default — Functionbat_default(f::Base.Callable, argname::Symbol, objectives...)
bat_default(f::Base.Callable, argname::Val, objectives...)Get the default value for argument argname of function f to use for objective(s).
objective(s) are mandatory arguments of function f that semantically constitute it's main objective(s), and that that a good default choice of optional arguments (e.g. choice of algorithm(s), etc.) may depend on. Which arguments are considered to be objectives is function-specific.
For example:
bat_default(bat_sample, :algorithm, density::PosteriorMeasure) == RandomWalk()
bat_default(bat_sample, Val(:algorithm), samples::DensitySampleVector) == OrderedResampling()BAT.bat_eff_sample_size — Functionbat_eff_sample_size(
v::Union{AbstractVector{<:Real},AbstractVectorOfSimilarVectors{<:Real}},
[algorithm::EffSampleSizeAlgorithm],
[context::BATContext]
)
bat_eff_sample_size(
smpls::DensitySampleVector,
[algorithm::EffSampleSizeAlgorithm],
[context::BATContext]
)Estimate effective sample size estimation for variate series v, resp. density samples smpls, separately for each degree of freedom.
Returns a NamedTuple of the shape
(result = eff_sample_size, ...)Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_findmedian — Functionbat_findmedian(
samples::DensitySampleVector
)The function computes the median of marginalized samples.
Returns a NamedTuple of the shape
(result = v, ...)Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_findmode — Functionbat_findmode(
target::BAT.AnySampleable,
[algorithm::BAT.AbstractModeEstimator],
[context::BATContext]
)::DensitySampleVectorEstimate the global mode of target.
Returns a NamedTuple of the shape
(result = X::DensitySampleVector, evaluated::EvaluatedMeasure, ...)(The field evaluated is only present if target is a measure.)
Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_initval — Functionbat_initval(
target::BAT.MeasureLike,
[algorithm::BAT.InitvalAlgorithm],
[context::BATContext]
)::V
bat_initval(
target::BAT.MeasureLike,
n::Integer,
[algorithm::BAT.InitvalAlgorithm],
[context::BATContext]
)::AbstractVector{<:V}Generate one or n random initial/starting value(s) suitable for target.
Assuming the variates of target are of type T, returns a NamedTuple of the shape
(result = X::AbstractVector{T}, ...)Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_integrate — Functionbat_integrate(
target::AnySampleable,
[algorithm::IntegrationAlgorithm],
[context::BATContext]
)::DensitySampleVectorCalculate the integral (evidence) of target.
Returns a NamedTuple of the shape
(result = X::Measurements.Measurement, evaluated::EvaluatedMeasure, ...)(The field evaluated is only present if target is a measure.)
Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_read — Functionbat_read(
filename::AbstractString,
[key,]
[algorithm::BATIOAlgorithm]
)Read data (optionally selected by key) from filename using algorithm.
Example:
smpls = bat_read("samples.hdf5", smpls).resultReturns (result = content, ...)
Result properties not listed here are specific to the output algorithm and are not part of the stable public API.
See bat_write.
Currently supported file formats are:
- HDF5 with file extension ".h5" or ".hdf5"
HDF5 I/O functionality is only available when the HDF5 package is loaded (e.g. via import HDF5).
BAT.bat_sample — Functionbat_sample(
target::BAT.AnySampleable,
[algorithm::BAT.AbstractSamplingAlgorithm],
[context::BATContext]
)::DensitySampleVectorDraw samples from target using algorithm.
Depending on sampling algorithm, the samples may be independent or correlated (e.g. when using MCMC).
Returns a NamedTuple of the shape
(result = X::DensitySampleVector, evaluated::EvaluatedMeasure, ...)(The field evaluated is only present if target is a measure.)
Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.bat_write — Functionbat_write(
filename::AbstractString,
content,
[algorithm::BATIOAlgorithm]
)Write content to file filename using algorithm.
Example:
smpls = bat_sample(posterior, ...).result
bat_write("samples.hdf5", smpls)Returns (result = filename, ...)
Result properties not listed here are specific to the output algorithm and are not part of the stable public API.
See bat_read.
Currently supported file formats are:
- HDF5 with file extension ".h5" or ".hdf5"
HDF5 I/O functionality is only available when the HDF5 package is loaded (e.g. via import HDF5).
BAT.bat_transform — Functionbat_transform(
how::AbstractTransformTarget,
object,
[algorithm::TransformAlgorithm]
)
bat_transform(
f,
object,
[algorithm::TransformAlgorithm]
)Transform object to another variate space defined/implied by target, res. using the transformation function f.
Returns a NamedTuple of the shape
(result = newdensity, f_transform = vartrafo::Function, ...)Result properties not listed here are algorithm-specific and are not part of the stable public API.
BAT.get_batcontext — Functionget_batcontext()::BATContext
get_batcontext(obj)::BATContextGets resp. sets the default computational context for BAT.
Will create and set a new default context if none exists.
Note: get_batcontext() does not have a stable return type. Code that needs type stability should pass a context to algorithms explicitly. BAT algorithms that call other algorithms must forward their context automatically, so context is always type stable within nested BAT algorithms.
See BATContext and set_batcontext.
BAT.set_batcontext — Functionset_batcontext(new_context::BATContext)
set_batcontext(;
precision = ...,
rng = ...,
ad = ...
)Sets the default computational context for BAT.
See BATContext and get_batcontext.
BAT.log_batdebug — Functionlog_batdebug(enable::Bool = true)Enable/disable debug-level logging for BAT and all BAT package extensions.
BAT.distbind — Functiondistbind(f_k, dist, ::typeof(merge))Performs a generalized monadic bind, in the functional programming sense, with a transition kernel f_k, a distribution dist, using merge to control the type of "flattening".
BAT.distprod — Functiondistprod(;a = some_dist, b = some_other_dist, ...)
distprod(();a = some_dist, b = some_other_dist, ...))
distprod([dist1, dist2, dist2, ...])Generate a product of distributions, returning either a distribution that has NamedTuples as variates, or arrays as variates.
BAT.lbqintegral — Functionlbqintegral(integrand, measure)
lbqintegral(likelihood, prior)Returns an object that represents the Lebesgue integral over a function in respect to s reference measure. It is also the non-normalized posterior measure that results from integrating the likelihood of a given observation in respect to a prior measure.
BAT.AbstractMCMCWeightingScheme — Typeabstract type AbstractMCMCWeightingScheme{T<:Real}Abstract class for weighting schemes for MCMC samples.
Weight values will have type T.
BAT.AbstractPosteriorMeasure — Typeabstract type AbstractPosteriorMeasure <: BATMeasure endAbstract type for posterior probability densities.
BAT.AbstractTransformTarget — Typeabstract type AbstractTransformTargetAbstract type for probability density transformation targets.
BAT.AdaptiveAffineTuning — Typestruct AdaptiveAffineTuning <: MCMCTransformTuningAdaptive cycle-based MCMC tuning strategy.
Adapts an affine space transformation based on the acceptance ratio and covariance of the previous samples.
Constructors:
AdaptiveAffineTuning(; fields...)
Fields:
λ::Float64: Controls the weight given to new covariance information in adapting the affine transform. Default: 0.5α::IntervalSets.ClosedInterval{Float64}: Metropolis-Hastings acceptance ratio target, tuning will try to adapt the affine transform to bring the acceptance ratio inside this interval. Default: ClosedInterval(0.15, 0.35)β::Float64: Controls how much the scale of the affine transform is widened/narrowed depending on the current MH acceptance ratio. Default: 1.5c::IntervalSets.ClosedInterval{Float64}: Interval for allowed scale of the affine transform distribution. Default: ClosedInterval(0.0001, 100.0)r::Real: Reweighting factor. Take accumulated sample statistics of previous tuning cycles into account with a relative weight ofr. Set to0to completely reset sample statistics between each tuning cycle. Default: 0.5
BAT.AssumeConvergence — Typestruct AssumeConvergence <: ConvergenceTest
No-op convergence algorithm for bat_convergence, will always declare convergence.
Constructors:
AssumeConvergence(converged::Bool = true)
Fields:
converged::Bool: Default: true
BAT.AutocorLenAlgorithm — Typeabstract type AutocorLenAlgorithmAbstract type for integrated autocorrelation length estimation algorithms.
BAT.BATContext — Typestruct BATContext{T}Experimental feature, not yet part of stable public API.
Set the default computational context for BAT.
Constructors:
BATContext{T}(rng::AbstractRNG, cunit::AbstractComputeUnit, ADSelector::AD)
BATContext(;
precision::Type{<:AbstractFloat} = ...,
rng::AbstractRNG = ...,
cunit::HeterogeneousComputing.AbstractComputeUnit = ...,
ad::Union{AutoDiffOperators.ADSelector, Module, Symbol, Val} = ...,
)See get_batcontext and set_batcontext.
BAT.BATHDF5IO — Typestruct BATHDF5IO <: BATIOAlgorithmSelects the BAT HDF5 format as the output format.
Constructors:
BATHDF5IO()
BAT.BATIOAlgorithm — Typeabstract type BATIOAlgorithmAbstract type for density transformation algorithms.
BAT.BinningAlgorithm — Typeabstract type BinningAlgorithmAbstract type for binning algorithms.
BAT.BrooksGelmanConvergence — Typestruct BrooksGelmanConvergence <: ConvergenceTestBrooks-Gelman maximum R^2 convergence test.
Constructors:
BrooksGelmanConvergence(; fields...)
Fields:
threshold::Float64: Default: 1.1corrected::Bool: Default: false
BAT.CuhreIntegration — Typestruct CuhreIntegration <: IntegrationAlgorithmCuhreIntegration integration algorithm.
Constructors:
CuhreIntegration(; fields...)
Fields:
pretransform::AbstractTransformTarget: Default: PriorToUniform()rtol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:RTOL))atol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:ATOL))minevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MINEVALS))maxevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXEVALS))key::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:KEY))nthreads::Int64: Default: Base.Threads.nthreads()strict::Bool: Default: true
This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).
BAT.DensitySample — Typestruct DensitySampleA weighted sample drawn according to an statistical density, e.g. a BAT.MeasureLike.
Constructors:
DensitySampleVector(v::Any, logd::Real, weight::Real, info::Any, aux::Any)
Fields:
v::Any: variate valuelogd::Real: log(density) value atvweight::Real: Weight of the sampleinfo::Any: Additional info on the provenance of the sample. Content depends on the sampling algorithm.aux::Any: Custom user-defined information attached to the sample.
Use DensitySampleVector to store vectors of multiple samples with an efficient column-based memory layout.
BAT.DensitySampleVector — Typestruct DensitySampleVector <: AbstractVector{<:DensitySample}A vector of DensitySample elements.
DensitySampleVector is currently a type alias for StructArrays.StructArray{<:DensitySample,...}, though this is subject to change without deprecation.
Constructors:
DensitySampleVector(
(
v::AbstractVector{<:AbstractVector{<:Real}},
logd::AbstractVector{<:Real},
weight::AbstractVector{<:Real},
info::AbstractVector{<:Any},
aux::AbstractVector{<:Any}
)
) DensitySampleVector(
v::AbstractVector,
logval::AbstractVector{<:Real};
weight::Union{AbstractVector{<:Real}, Symbol} = fill(1, length(eachindex(v))),
info::AbstractVector = fill(nothing, length(eachindex(v))),
aux::AbstractVector = fill(nothing, length(eachindex(v)))
)With weight = :multiplicity repeated samples will be replaced by a single sample, with a weight equal to the number of repetitions.
BAT.DivonneIntegration — Typestruct DivonneIntegration <: IntegrationAlgorithmDivonneIntegration integration algorithm.
Constructors:
DivonneIntegration(; fields...)
Fields:
pretransform::AbstractTransformTarget: Default: PriorToUniform()rtol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:RTOL))atol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:ATOL))minevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MINEVALS))maxevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXEVALS))key1::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:KEY1))key2::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:KEY2))key3::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:KEY3))maxpass::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXPASS))border::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:BORDER))maxchisq::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXCHISQ))mindeviation::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MINDEVIATION))ngiven::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NGIVEN))ldxgiven::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:LDXGIVEN))nextra::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NEXTRA))nthreads::Int64: Default: Base.Threads.nthreads()strict::Bool: Default: true
This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).
BAT.DoNotTransform — Typestruct DoNotTransform <: AbstractTransformTargetThe identity density transformation target, specifies that densities should not be transformed.
Constructors:
DoNotTransform()
BAT.EffSampleSizeAlgorithm — Typeabstract type EffSampleSizeAlgorithmAbstract type for integrated autocorrelation length estimation algorithms.
BAT.EffSampleSizeFromAC — Typestruct EffSampleSizeFromAC <: EffSampleSizeAlgorithmEffective sample size estimation based on the integrated autocorrelation length of the samples.
Constructors:
EffSampleSizeFromAC(; fields...)
Fields:
acalg::AutocorLenAlgorithm: Default: GeyerAutocorLen()
BAT.EvaluatedMeasure — Typestruct EvaluatedMeasure <: BATMeasureCombined a measure with samples, and other information on it.
Constructors:
em = EvaluatedMeasure(
measure;
samples = ..., approx = ..., mass = ..., mode = ...,
_generator = ...
)
BAT.unevaluated(em) === measureBAT.ExplicitInit — Typestruct ExplicitInit <: InitvalAlgorithmUses initial values from a given vector of one or more values/variates. The values are used in the order they appear in the vector, not randomly.
Constructors:
ExplicitInit(; fields...)
Fields:
xs::AbstractVector
BAT.FixedMGVISchedule — Typeabstract type FixedMGVISchedule <: BAT.MGVIScheduleAbstract supertype for MGVI sampling schedules.
Constructors:
FixedMGVISchedule(; fields...)
Fields:
nsamples::AbstractVector{<:Real}: Default: range(12, 1000, length = 10)
Constructors:
FixedMGVISchedule(nsamples::AbstractVector{<:Real}): The number of samples to draw at each MGVI step. The length ofnsamplesimplies the total number of steps. The number of samples will be rounded to integer values if necessary, to allow for constructions likeFixedMGVISchedule(range(12, 1000, length = 10)).
Fields:
nsamples::AbstractVector{<:Real}: See constructor above.
See MGVISampling.
BAT.FixedNBins — TypeFixedNBins(nbins::Int)Selects a fixed number of bins.
Constructor: FixedNBins(; fields...)
Fields:
nbins::Int64: Default: 200
BAT.FreedmanDiaconisBinning — Typestruct FreedmanDiaconisBinning <: BinningAlgorithmSelects automatic binning based on the Freedman–Diaconis rule.
Constructor: FreedmanDiaconisBinning()
BAT.GelmanRubinConvergence — Typestruct GelmanRubinConvergence <: ConvergenceTestGelman-Rubin maximum R^2 convergence test.
Constructors:
GelmanRubinConvergence(; fields...)
Fields:
threshold::Float64: Default: 1.1
BAT.GeyerAutocorLen — Typestruct GeyerAutocorLen <: AutocorLenAlgorithmIntegrated autocorrelation length estimation based on Geyer’s initial monotone sequence criterion
See C. J. Geyer, "Practical Markov Chain Monte Carlo" (1992) and C. J. Geyer, "Introduction to Markov Chain Monte Carlo" (2011).
Constructors:
GeyerAutocorLen()
The same algorithm is used by STAN (v2.21) and MCMCChains.jl (v3.0, function ess_rhat).
BAT.HamiltonianMC — Typestruct HamiltonianMC <: MCMCAlgorithmThe Hamiltonian Monte Carlo (HMC) sampling algorithm.
Uses the HMC implementation provided by the package AdvancedHMC.
HMC uses gradients of the target measure's density, so your BATContext needs to include an ADSelector to specify which automatic differentiation backend should be used.
- Note: The fields of
HamiltonianMCare still subject to change, and not
yet part of stable public BAT API!*
Constructors:
HamiltonianMC(; fields...)
Fields:
metric::BAT.HMCMetric: Default: UnitEuclideanMetric()integrator::Any: Default: extdefault(pkgext(Val(:AdvancedHMC)), Val(:DEFAULTINTEGRATOR))termination::Any: Default: extdefault(pkgext(Val(:AdvancedHMC)), Val(:DEFAULTTERMINATION_CRITERION))
BAT.IdentityTransformAlgorithm — Typestruct IdentityTransformAlgorithm <: TransformAlgorithmA no-op density transform algorithm that leaves any density unchanged.
Constructors:
IdentityTransformAlgorithm()
BAT.IIDSampling — Typestruct IIDSampling <: AbstractSamplingAlgorithmSample via Random.rand.
Constructors:
IIDSampling(; fields...)
Fields:
nsamples::Int64: Default: 10 ^ 5
BAT.InitFromIID — Typestruct InitFromIID <: InitvalAlgorithmGenerates initial values for sampling, optimization, etc. by random resampling from a given set of samples.
Constructors:
InitFromIID()
BAT.InitFromSamples — Typestruct InitFromSamples <: InitvalAlgorithmGenerates initial values for sampling, optimization, etc. by direct sampling from a given i.i.d. sampleable source.
Constructors:
InitFromSamples()
BAT.InitFromTarget — Typestruct InitFromTarget <: InitvalAlgorithmGenerates initial values for sampling, optimization, etc. by direct i.i.d. sampling a suitable component of that target density (e.g. it's prior) that supports it.
If the target is supports direct i.i.d. sampling, e.g. because it is a distribution, initial values are sampled directly from the target.
If the target is a posterior density, initial values are sampled from the prior (or the prior's prior if the prior is a posterior itself, etc.).
If the target is a sampled density, initial values are (re-)sampled from the available samples.
Constructors:
InitFromTarget()
BAT.InitvalAlgorithm — Typeabstract type BAT.InitvalAlgorithmAbstract type for BAT initial/starting value generation algorithms.
Many algorithms in BAT, like MCMC and optimization, need initial/starting values.
BAT.IntegrationAlgorithm — Typeabstract type IntegrationAlgorithmAbstract type for integration algorithms.
BAT.KishESS — Typestruct KishESS <: EffSampleSizeAlgorithmKish's effective sample size estimator, uses only the sample weights.
Constructors:
KishESS()
BAT.MaxDensitySearch — TypeMaxDensitySearch <: AbstractModeEstimatorConstructors:
MaxDensitySearch()Estimate the mode as the variate with the highest posterior density value within a given set of samples.
BAT.MCMCAlgorithm — Typeabstract type MCMCAlgorithmAbstract type for Markov chain Monte Carlo algorithms.
To implement a new MCMC algorithm, subtypes of both MCMCAlgorithm and MCMCChainState are required.
BAT.MCMCBurninAlgorithm — Typeabstract type MCMCBurninAlgorithmAbstract type for MCMC burn-in algorithms.
BAT.MCMCChainPoolInit — Typestruct MCMCChainPoolInit <: MCMCInitAlgorithmMCMC chain pool initialization strategy.
Constructors:
MCMCChainPoolInit(; fields...)
Fields:
init_tries_per_chain::IntervalSets.ClosedInterval{Int64}: Default: ClosedInterval(8, 128)nsteps_init::Int64: Default: 1000initval_alg::InitvalAlgorithm: Default: InitFromTarget()strict::Bool: Default: true
BAT.MCMCRetryInit — Typestruct MCMCRetryInit <: MCMCInitAlgorithmTODO
Constructors:
MCMCRetryInit(; fields...)
Fields:
max_init_tries::Int64: Default: 20nsteps_init::Int64: Default: 250initval_alg::InitvalAlgorithm: Default: InitFromTarget()strict::Bool: Default: true
BAT.MCMCInitAlgorithm — Typeabstract type MCMCInitAlgorithmAbstract type for MCMC initialization algorithms.
BAT.MCMCMultiCycleBurnin — Typestruct MCMCMultiCycleBurnin <: MCMCBurninAlgorithmA multi-cycle MCMC burn-in algorithm.
Constructors:
MCMCMultiCycleBurnin(; fields...)
Fields:
nsteps_per_cycle::Int64: Default: 10000max_ncycles::Int64: Default: 30nsteps_final::Int64: Default: div(nstepspercycle, 10)
BAT.MCMCProposalTuning — Typeabstract type MCMCProposalTuningAbstract type for MCMC tuning algorithms.
BAT.MCMCTransformTuning — Typeabstract type MCMCTransformTuningAbstract type for MCMC tuning algorithms.
BAT.MGVISampling — Typestruct MGVISampling <: AbstractUltraNestAlgorithmReactivSamples via Metric Gaussian Variational Inference, using the MGVI.jl Julia implementation of the algorithm.
Constructors:
MGVISampling(; fields...)
Fields:
pretransform::AbstractTransformTarget: Pre-transformation to apply to the target measure before sampling.nsamples::Int: Number is independent samples to draw. MGVI will generate symmetical samples, so it will generate2*nsamplessamples in total, but onlynsamples` independent samples.schedule::MGVISchedule: MGVI schedule, by default aFixedMGVISchedule.config::MGVI.MGVIConfig: MGVI configuration.
This functionality is only available when the package MGVI is loaded (e.g. via import MGVI).
BAT.ModeAsDefined — Typestruct ModeAsDefined <: AbstractModeEstimatorGet the mode as defined by the density, resp. the underlying distribution (if available), via StatsBase.mode.
Constructors:
ModeAsDefined()
BAT.NoMCMCProposalTuning — TypeNoMCMCProposalTuning <: MCMCProposalTuningDo not perform any MCMC proposal tuning.
BAT.NoMCMCTransformTuning — TypeNoMCMCTransformTuning <: MCMCTransformTuningDo not perform any MCMC transform turing.
BAT.OptimAlg — TypeOptimAlgSelects an optimization algorithm from the Optim.jl package.
Note that when using first order algorithms like Optim.LBFGS, your BATContext needs to include an ADSelector that specifies which automatic differentiation backend should be used.
Constructors:
OptimAlg(; fields...)
optimalg must be an Optim.AbstractOptimizer.
Fields:
optalg::Any: Default: extdefault(pkgext(Val(:Optim)), Val(:DEFAULTOPTALG))pretransform::AbstractTransformTarget: Default: PriorToNormal()init::InitvalAlgorithm: Default: InitFromTarget()maxiters::Int64: Default: 1000maxtime::Float64: Default: NaNabstol::Float64: Default: NaNreltol::Float64: Default: 0.0kwargs::NamedTuple: Default: (;)
BAT.OptimizationAlg — Typestruct OptimizationAlgSelects an optimization algorithm from the Optimization.jl package. Note that when using first order algorithms like OptimizationOptimJL.LBFGS, your BATContext needs to have ad set to an automatic differentiation backend.
Constructors:
OptimizationAlg(; fields...)
optalg must be an Optimization.AbstractOptimizer. The field kwargs can be used to pass additional keywords to the optimizers See the Optimization.jl documentation for the available keyword arguments. Fields:
optalg::Any: Default: extdefault(pkgext(Val(:Optimization)), Val(:DEFAULTOPTALG))pretransform::AbstractTransformTarget: Default: PriorToNormal()init::InitvalAlgorithm: Default: InitFromTarget()maxiters::Int64: Default: 1000maxtime::Float64: Default: NaNabstol::Float64: Default: NaNreltol::Float64: Default: 0.0kwargs::NamedTuple: Default: (;)
BAT.OrderedResampling — Typestruct OrderedResampling <: AbstractSamplingAlgorithmEfficiently resamples from a given series of samples, keeping the order of samples.
Can be used to efficiently convert weighted samples into samples with unity weights.
Constructors:
OrderedResampling(; fields...)
Fields:
nsamples::Int64: Default: 10 ^ 5
BAT.PosteriorMeasure — Typestruct PosteriorMeasure{Li,Pr<:AbstractMeasure} <: AbstractPosteriorMeasureA representation of a PosteriorMeasure, based a likelihood and prior. Likelihood and prior be accessed via
getlikelihood(posterior::PosteriorMeasure)::Li
getprior(posterior::PosteriorMeasure)::PrConstructors:
PosteriorMeasure(likelihood, prior)PosteriorMeasure{T<:Real}(likelihood, prior)
Fields:
likelihood::Anyprior::MeasureBase.AbstractMeasure
BAT.PriorSubstitution — Typestruct PriorSubstitution <: TransformAlgorithmSubstitute the prior by a given distribution and transform the likelihood accordingly. The log(abs(jacobian)) of the transformation does not need to be auto-differentiable even for operations that use the gradient of the posterior.
Constructors:
PriorSubstitution()
BAT.PriorToNormal — Typestruct PriorToNormal <: AbstractTransformTargetSpecifies that posterior densities should be transformed in a way that makes their pior equivalent to a standard multivariate normal distribution with an identity covariance matrix.
Constructors:
PriorToNormal()
BAT.PriorToUniform — Typestruct PriorToUniform <: AbstractTransformTargetSpecifies that posterior densities should be transformed in a way that makes their pior equivalent to a uniform distribution over the unit hypercube.
Constructors:
PriorToUniform()
BAT.RAMTuning — Typestruct RAMTuning <: MCMCTransformTuningTunes MCMC spaces transformations based on the Robust adaptive Metropolis algorithm.
In constrast to the original RAM algorithm, RAMTuning does not use the covariance estimate to change a proposal distribution, but instead uses it as the bases for an affine transformation. The sampling process is mathematically equivalent, though.
Constructors:
RAMTuning(; fields...)
Fields:
target_acceptance::Float64: MCMC target acceptance ratio. Default: 0.234σ_target_acceptance::Float64: Width aroundtarget_acceptance. Default: 0.05gamma::Float64: Negative adaption rate exponent. Default: 2 / 3
BAT.RandomWalk — Typestruct RandomWalk <: MCMCAlgorithmMetropolis-Hastings MCMC sampling algorithm.
Constructors:
RandomWalk(; fields...)
Fields:
proposaldist::Union{MeasureBase.AbstractMeasure, Distributions.ContinuousDistribution{<:Union{Distributions.Univariate, Distributions.Multivariate}}}: Default: TDist(1.0)
BAT.RandResampling — Typestruct RandResampling <: AbstractSamplingAlgorithmResamples from a given set of samples.
Constructors:
RandResampling(; fields...)
Fields:
nsamples::Int64: Default: 10 ^ 5
BAT.RepetitionWeighting — Typestruct RepetitionWeighting{T<:AbstractFloat} <: AbstractMCMCWeightingScheme{T}Sample weighting scheme suitable for sampling algorithms which may repeated samples multiple times in direct succession (e.g. RandomWalk). The repeated sample is stored only once, with a weight equal to the number of times it has been repeated (e.g. because a Markov chain has not moved during a sampling step).
Constructors:
RepetitionWeighting()
BAT.RiceBinning — Typestruct RiceBinning <: BinningAlgorithmSelects automatic binning based on the Rice rule.
Constructor: RiceBinning()
BAT.SampleMedianEstimator — Typestruct SampleMedianEstimator <: AbstractMedianEstimatorGet median values from samples using standard Julia statistics functions.
Constructors:
SampleMedianEstimator()
BAT.ScottBinning — Typestruct ScottBinning <: BinningAlgorithmSelects automatic binning based on Scott's normal reference rule.
Constructor: ScottBinning()
BAT.SokalAutocorLen — Typestruct SokalAutocorLen <: AutocorLenAlgorithmIntegrated autocorrelation length estimation based on the automated windowing procedure descibed in A. D. Sokal, "Monte Carlo Methods in Statistical Mechanics" (1996)
Same procedure is used by the emcee Python package (v3.0).
Constructors:
SokalAutocorLen(; fields...)
Fields:
c::Int64: Step size for window search Default: 5
BAT.SquareRootBinning — Typestruct SquareRootBinning <: BinningAlgorithmSelects automatic binning based on the Square-root choice.
Constructor: SquareRootBinning()
BAT.SturgesBinning — Typestruct SturgesBinning <: BinningAlgorithmSelects automatic binning based on Sturges' formula.
Constructor: SturgesBinning()
BAT.SuaveIntegration — Typestruct SuaveIntegration <: IntegrationAlgorithmSuaveIntegration integration algorithm.
Constructors:
SuaveIntegration(; fields...)
Fields:
pretransform::AbstractTransformTarget: Default: PriorToUniform()rtol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:RTOL))atol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:ATOL))minevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MINEVALS))maxevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXEVALS))nnew::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NNEW))nmin::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NMIN))flatness::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:FLATNESS))nthreads::Int64: Default: Base.Threads.nthreads()strict::Bool: Default: true
This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).
BAT.ToRealVector — Typestruct ToRealVector <: AbstractTransformTargetSpecifies that the input should be transformed into a measure over the space of real-valued flat vectors.
Constructors:
ToRealVector()
BAT.TransformAlgorithm — Typeabstract type TransformAlgorithmAbstract type for density transformation algorithms.
BAT.TransformedMCMC — Typestruct TransformedMCMC <: AbstractSamplingAlgorithmSamples a probability density using Markov chain Monte Carlo.
Constructors:
TransformedMCMC(; fields...)
Fields:
proposal::BAT.MCMCProposal: Default: RandomWalk(proposaldist = TDist(1.0))proposal_tuning::MCMCProposalTuning: Default: batdefault(TransformedMCMC, Val(:proposaltuning), proposal)pretransform::AbstractTransformTarget: Default: bat_default(TransformedMCMC, Val(:pretransform), proposal)adaptive_transform::BAT.AbstractAdaptiveTransform: Default: batdefault(TransformedMCMC, Val(:adaptivetransform), proposal)transform_tuning::MCMCTransformTuning: Default: batdefault(TransformedMCMC, Val(:transformtuning), adaptive_transform)tempering::MCMCTempering: Default: bat_default(TransformedMCMC, Val(:tempering), proposal)nchains::Int64: Default: 4nwalkers::Int64: Default: batdefault(TransformedMCMC, Val(:nwalkers), proposal, pretransform, transformtuning, nchains)nsteps::Int64: Default: batdefault(TransformedMCMC, Val(:nsteps), proposal, pretransform, transformtuning, nchains, nwalkers)init::MCMCInitAlgorithm: Default: batdefault(TransformedMCMC, Val(:init), proposal, pretransform, transformtuning, nchains, nwalkers, nsteps)burnin::MCMCBurninAlgorithm: Default: batdefault(TransformedMCMC, Val(:burnin), proposal, pretransform, transformtuning, nchains, nwalkers, nsteps)convergence::BAT.ConvergenceTest: Default: BrooksGelmanConvergence()strict::Bool: Default: truestore_burnin::Bool: Default: falsenonzero_weights::Bool: Default: truesample_weighting::AbstractMCMCWeightingScheme: Default: RepetitionWeighting()callback::Function: Default: nop_func
BAT.VEGASIntegration — Typestruct VEGASIntegration <: IntegrationAlgorithmVEGASIntegration integration algorithm.
Constructors:
VEGASIntegration(; fields...)
Fields:
pretransform::AbstractTransformTarget: Default: PriorToUniform()rtol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:RTOL))atol::Float64: Default: ext_default(pkgext(Val(:Cuba)), Val(:ATOL))minevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MINEVALS))maxevals::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:MAXEVALS))nstart::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NSTART))nincrease::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NINCREASE))nbatch::Int64: Default: ext_default(pkgext(Val(:Cuba)), Val(:NBATCH))nthreads::Int64: Default: Base.Threads.nthreads()strict::Bool: Default: true
This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).
BAT.unevaluated — FunctionBAT.unevaluated(obj)If obj is an evaluated object, like a EvaluatedMeasure, return the original (unevaluated) object. Otherwise, return obj.
BAT.AbstractMedianEstimator — Typeabstract type BAT.AbstractMedianEstimatorAbstract type for BAT optimization algorithms.
A typical application for optimization in BAT is mode estimation (see bat_findmode),
BAT.AbstractModeEstimator — Typeabstract type BAT.AbstractModeEstimatorAbstract type for BAT optimization algorithms.
A typical application for optimization in BAT is mode estimation (see bat_findmode),
BAT.AbstractSamplingAlgorithm — Typeabstract type BAT.AbstractSamplingAlgorithmAbstract type for BAT sampling algorithms. See bat_sample.
BAT.ConvergenceTest — Typeabstract type ConvergenceTestAbstract type for integrated autocorrelation length estimation algorithms.
BAT.MGVISchedule — Type