API Documentation

This is the stable public API of BAT. Forward/backward compatibility follows Julia's semantic versioning rules.

Types

Functions and macros

Documentation

BAT.bat_convergenceFunction
bat_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.

Note

Do not add add methods to bat_convergence, add methods to bat_convergence_impl instead.

source
BAT.bat_defaultFunction
bat_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) == MetropolisHastings()
bat_default(bat_sample, Val(:algorithm), samples::DensitySampleVector) == OrderedResampling()
source
BAT.bat_eff_sample_sizeFunction
bat_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.

Note

Do not add add methods to bat_eff_sample_size, add methods to bat_eff_sample_size_impl instead.

source
BAT.bat_findmedianFunction
bat_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.

Note

Do not add add methods to bat_findmedian, add methods to bat_findmedian_impl instead.

source
BAT.bat_findmodeFunction
bat_findmode(
    target::BAT.AnySampleable,
    [algorithm::BAT.AbstractModeEstimator],
    [context::BATContext]
)::DensitySampleVector

Estimate the global mode of target.

Returns a NamedTuple of the shape

(result = X::DensitySampleVector, ...)

Result properties not listed here are algorithm-specific and are not part of the stable public API.

Note

Do not add add methods to bat_findmode, add methods to bat_findmode_impl instead.

source
BAT.bat_initvalFunction
bat_initval(
    target::BAT.AnyMeasureOrDensity,
    [algorithm::BAT.InitvalAlgorithm],
    [context::BATContext]
)::V

bat_initval(
    target::BAT.AnyMeasureOrDensity,
    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.

Note

Do not add add methods to bat_initval, add methods like

bat_initval_impl(target::AnyMeasureOrDensity, algorithm::InitvalAlgorithm, context::BATContext)
bat_initval_impl(target::AnyMeasureOrDensity, n::Integer, algorithm::InitvalAlgorithm, context::BATContext)

to bat_initval_impl instead.

source
BAT.bat_integrateFunction
bat_integrate(
    target::AnySampleable,
    [algorithm::IntegrationAlgorithm],
    [context::BATContext]
)::DensitySampleVector

Calculate the integral (evidence) of target.

Returns a NamedTuple of the shape

(result = X::Measurements.Measurement, ...)

Result properties not listed here are algorithm-specific and are not part of the stable public API.

Note

Do not add add methods to bat_integrate, add methods to bat_integrate_impl instead.

source
BAT.bat_readFunction
bat_read(
    filename::AbstractString,
    [key,]
    [algorithm::BATIOAlgorithm]
)::AbstractMeasureOrDensity

Read data (optionally selected by key) from filename using algorithm.

Example:

smpls = bat_read("samples.hdf5", smpls).result

Returns (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"
Note

HDF5 I/O functionality is only available when the HDF5 package is loaded (e.g. via import HDF5).

Note

Do not add add algorithms to bat_read, add algorithms to bat_read_impl instead.

source
BAT.bat_sampleFunction
bat_sample(
    target::BAT.AnySampleable,
    [algorithm::BAT.AbstractSamplingAlgorithm],
    [context::BATContext]
)::DensitySampleVector

Draw 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, ...)

Result properties not listed here are algorithm-specific and are not part of the stable public API.

Note

Do not add add methods to bat_sample, add methods to bat_sample_impl instead.

source
BAT.bat_writeFunction
bat_write(
    filename::AbstractString,
    content,
    [algorithm::BATIOAlgorithm]
)::AbstractMeasureOrDensity

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"
Note

HDF5 I/O functionality is only available when the HDF5 package is loaded (e.g. via import HDF5).

Note

Do not add add algorithms to bat_write, add algorithms to bat_write_impl instead.

source
BAT.bat_transformFunction
bat_transform(
    how::AbstractTransformTarget,
    object,
    [algorithm::TransformAlgorithm]
)::AbstractMeasureOrDensity

bat_transform(
    f,
    object,
    [algorithm::TransformAlgorithm]
)::AbstractMeasureOrDensity

Transform object to another variate space defined/implied by target, res. using the transformation function f.

Returns a NamedTuple of the shape

(result = newdensity::AbstractMeasureOrDensity, trafo = vartrafo::Function, ...)

Result properties not listed here are algorithm-specific and are not part of the stable public API.

Note

As a convenience,

flat_smpls, f_flatten = bat_transform(Vector, measure)
flat_smpls, f_flatten = bat_transform(Vector, samples)

can be used to flatten a the variate type of a measure (res. samples of a measure) to something like Vector{<:Real}.

source
BAT.get_batcontextFunction
get_batcontext()::BATContext
get_batcontext(obj)::BATContext

Gets 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.

source
BAT.distbindFunction
distbind(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".

source
BAT.distprodFunction
distprod(;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.

source
BAT.lbqintegralFunction
lbqintegral(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.

source
BAT.AdaptiveMHTuningType
struct AdaptiveMHTuning <: MHProposalDistTuning

Adaptive MCMC tuning strategy for Metropolis-Hastings samplers.

Adapts the proposal function based on the acceptance ratio and covariance of the previous samples.

Constructors:

  • AdaptiveMHTuning(; fields...)

Fields:

  • λ::Float64: Controls the weight given to new covariance information in adapting the proposal distribution. Default: 0.5

  • α::IntervalSets.ClosedInterval{Float64}: Metropolis-Hastings acceptance ratio target, tuning will try to adapt the proposal distribution to bring the acceptance ratio inside this interval. Default: ClosedInterval(0.15, 0.35)

  • β::Float64: Controls how much the spread of the proposal distribution is widened/narrowed depending on the current MH acceptance ratio. Default: 1.5

  • c::IntervalSets.ClosedInterval{Float64}: Interval for allowed scale/spread of the proposal 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 of r. Set to 0 to completely reset sample statistics between each tuning cycle. Default: 0.5

source
BAT.AssumeConvergenceType

struct AssumeConvergence <: ConvergenceTest

No-op convergence algorithm for bat_convergence, will always declare convergence.

Constructors:

  • AssumeConvergence(converged::Bool = true)

Fields:

  • converged::Bool: Default: true
source
BAT.AutocorLenAlgorithmType
abstract type AutocorLenAlgorithm

Abstract type for integrated autocorrelation length estimation algorithms.

source
BAT.BATContextType
struct 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::AutoDiffOperators.ADSelector = ...,
)

See get_batcontext and set_batcontext.

source
BAT.BrooksGelmanConvergenceType
struct BrooksGelmanConvergence <: ConvergenceTest

Brooks-Gelman maximum R^2 convergence test.

Constructors:

  • BrooksGelmanConvergence(; fields...)

Fields:

  • threshold::Float64: Default: 1.1

  • corrected::Bool: Default: false

source
BAT.CuhreIntegrationType
struct CuhreIntegration <: IntegrationAlgorithm

CuhreIntegration integration algorithm.

Constructors:

  • CuhreIntegration(; fields...)

Fields:

  • trafo::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

Note

This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).

source
BAT.DensitySampleType
struct DensitySample

A weighted sample drawn according to an statistical density, e.g. a BAT.AbstractMeasureOrDensity.

Constructors:

  • DensitySampleVector(v::Any, logd::Real, weight::Real, info::Any, aux::Any)

Fields:

  • v::Any: variate value

  • logd::Real: log(density) value at v

  • weight::Real: Weight of the sample

  • info::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.

source
BAT.DensitySampleVectorType
struct 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.

source
BAT.DivonneIntegrationType
struct DivonneIntegration <: IntegrationAlgorithm

DivonneIntegration integration algorithm.

Constructors:

  • DivonneIntegration(; fields...)

Fields:

  • trafo::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

Note

This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).

source
BAT.DoNotTransformType
struct DoNotTransform <: AbstractTransformTarget

The identity density transformation target, specifies that densities should not be transformed.

Constructors:

  • DoNotTransform()
source
BAT.EffSampleSizeFromACType
struct EffSampleSizeFromAC <: EffSampleSizeAlgorithm

Effective sample size estimation based on the integrated autocorrelation length of the samples.

Constructors:

  • EffSampleSizeFromAC(; fields...)

Fields:

  • acalg::AutocorLenAlgorithm: Default: GeyerAutocorLen()
source
BAT.EvaluatedMeasureType
struct EvaluatedMeasure <: BATMeasure

Combined a measure with samples, and other information on it.

Constructors:

EvaluatedMeasure(
    measure;
    samples = ..., approx = ..., mass = ..., mode = ...,
    _generator = ...
)
Note

Field _generator does not form part of the stable public API and is subject to change without deprecation.

source
BAT.ExplicitInitType
struct ExplicitInit <: InitvalAlgorithm

Uses 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
source
BAT.GelmanRubinConvergenceType
struct GelmanRubinConvergence <: ConvergenceTest

Gelman-Rubin maximum R^2 convergence test.

Constructors:

  • GelmanRubinConvergence(; fields...)

Fields:

  • threshold::Float64: Default: 1.1
source
BAT.HamiltonianMCType
struct HamiltonianMC <: MCMCAlgorithm

The 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 HamiltonianMC are still subject to change, and not

yet part of stable public BAT API!*

Constructors:

  • HamiltonianMC(; fields...)

Fields:

  • metric::BAT.HMCMetric: Default: DiagEuclideanMetric()

  • integrator::Any: Default: extdefault(pkgext(Val(:AdvancedHMC)), Val(:DEFAULTINTEGRATOR))

  • termination::Any: Default: extdefault(pkgext(Val(:AdvancedHMC)), Val(:DEFAULTTERMINATION_CRITERION))

  • tuning::BAT.HMCTuningAlgorithm: Default: StanHMCTuning()

Note

HamiltonianMC is only available if the AdvancedHMC package is loaded (e.g. via import AdvancedHMC).

source
BAT.IdentityTransformAlgorithmType
struct IdentityTransformAlgorithm <: TransformAlgorithm

A no-op density transform algorithm that leaves any density unchanged.

Constructors:

  • IdentityTransformAlgorithm()
source
BAT.IIDSamplingType
struct IIDSampling <: AbstractSamplingAlgorithm

Sample via Random.rand. Only supported for posteriors of type Distributions.MultivariateDistribution and BAT.DistLikeMeasure.

Constructors:

  • IIDSampling(; fields...)

Fields:

  • nsamples::Int64: Default: 10 ^ 5
source
BAT.InitFromIIDType
struct InitFromIID <: InitvalAlgorithm

Generates initial values for sampling, optimization, etc. by random resampling from a given set of samples.

Constructors:

  • InitFromIID()
source
BAT.InitFromSamplesType
struct InitFromSamples <: InitvalAlgorithm

Generates initial values for sampling, optimization, etc. by direct sampling from a given i.i.d. sampleable source.

Constructors:

  • InitFromSamples()
source
BAT.InitFromTargetType
struct InitFromTarget <: InitvalAlgorithm

Generates 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()
source
BAT.InitvalAlgorithmType
abstract type BAT.InitvalAlgorithm

Abstract type for BAT initial/starting value generation algorithms.

Many algorithms in BAT, like MCMC and optimization, need initial/starting values.

source
BAT.KishESSType
struct KishESS <: EffSampleSizeAlgorithm

Kish's effective sample size estimator, uses only the sample weights.

Constructors:

  • KishESS()
source
BAT.MaxDensitySearchType
MaxDensitySearch <: AbstractModeEstimator

Constructors:

MaxDensitySearch()

Estimate the mode as the variate with the highest posterior density value within a given set of samples.

source
BAT.MCMCAlgorithmType
abstract type MCMCAlgorithm

Abstract type for Markov chain Monte Carlo algorithms.

To implement a new MCMC algorithm, subtypes of both MCMCAlgorithm and MCMCIterator are required.

Note

The details of the MCMCIterator and MCMCAlgorithm API required to implement a new MCMC algorithm currently do not (yet) form part of the stable API and are subject to change without deprecation.

source
BAT.MCMCChainPoolInitType
struct MCMCChainPoolInit <: MCMCInitAlgorithm

MCMC chain pool initialization strategy.

Constructors:

  • MCMCChainPoolInit(; fields...)

Fields:

  • init_tries_per_chain::IntervalSets.ClosedInterval{Int64}: Default: ClosedInterval(8, 128)

  • nsteps_init::Int64: Default: 1000

  • initval_alg::InitvalAlgorithm: Default: InitFromTarget()

source
BAT.MCMCMultiCycleBurninType
struct MCMCMultiCycleBurnin <: MCMCBurninAlgorithm

A multi-cycle MCMC burn-in algorithm.

Constructors:

  • MCMCMultiCycleBurnin(; fields...)

Fields:

  • nsteps_per_cycle::Int64: Default: 10000

  • max_ncycles::Int64: Default: 30

  • nsteps_final::Int64: Default: div(nstepspercycle, 10)

source
BAT.MCMCNoOpTuningType
MCMCNoOpTuning <: MCMCTuningAlgorithm

No-op tuning, marks MCMC chains as tuned without performing any other changes on them. Useful if chains are pre-tuned or tuning is an internal part of the MCMC sampler implementation.

source
BAT.MCMCSamplingType
struct MCMCSampling <: AbstractSamplingAlgorithm

Samples a probability density using Markov chain Monte Carlo.

Constructors:

  • MCMCSampling(; fields...)

Fields:

  • mcalg::MCMCAlgorithm: Default: MetropolisHastings()

  • trafo::AbstractTransformTarget: Default: bat_default(MCMCSampling, Val(:trafo), mcalg)

  • nchains::Int64: Default: 4

  • nsteps::Int64: Default: bat_default(MCMCSampling, Val(:nsteps), mcalg, trafo, nchains)

  • init::MCMCInitAlgorithm: Default: bat_default(MCMCSampling, Val(:init), mcalg, trafo, nchains, nsteps)

  • burnin::MCMCBurninAlgorithm: Default: bat_default(MCMCSampling, Val(:burnin), mcalg, trafo, nchains, nsteps)

  • convergence::BAT.ConvergenceTest: Default: BrooksGelmanConvergence()

  • strict::Bool: Default: true

  • store_burnin::Bool: Default: false

  • nonzero_weights::Bool: Default: true

  • callback::Function: Default: nop_func

source
BAT.MetropolisHastingsType
struct MetropolisHastings <: MCMCAlgorithm

Metropolis-Hastings MCMC sampling algorithm.

Constructors:

  • MetropolisHastings(; fields...)

Fields:

  • proposal::BAT.ProposalDistSpec: Default: MvTDistProposal()

  • weighting::AbstractMCMCWeightingScheme: Default: RepetitionWeighting()

  • tuning::MHProposalDistTuning: Default: AdaptiveMHTuning()

source
BAT.MHProposalDistTuningType
abstract type MHProposalDistTuning

Abstract type for Metropolis-Hastings tuning strategies for proposal distributions.

source
BAT.ModeAsDefinedType
struct ModeAsDefined <: AbstractModeEstimator

Get the mode as defined by the density, resp. the underlying distribution (if available), via StatsBase.mode.

Constructors:

  • ModeAsDefined()
source
BAT.OptimAlgType
OptimAlg

Selects 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))

  • trafo::AbstractTransformTarget: Default: PriorToGaussian()

  • init::InitvalAlgorithm: Default: InitFromTarget()

Note

This algorithm is only available if the Optim package is loaded (e.g. via import Optim.

source
BAT.OrderedResamplingType
struct OrderedResampling <: AbstractSamplingAlgorithm

Efficiently 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
source
BAT.PosteriorMeasureType
struct PosteriorMeasure{
    Li<:AbstractMeasureOrDensity,
    Pr<:DistLikeMeasure,
    ...
} <: AbstractPosteriorMeasure

A representation of a PosteriorMeasure, based a likelihood and prior. Likelihood and prior be accessed via

getlikelihood(posterior::PosteriorMeasure)::Li
getprior(posterior::PosteriorMeasure)::Pr

Constructors:

  • PosteriorMeasure(likelihood, prior)
  • PosteriorMeasure{T<:Real}(likelihood, prior)

likelihood and prior must be convertible to an AbstractMeasureOrDensity.

Fields:

  • likelihood::BAT.AbstractMeasureOrDensity

  • prior::BAT.AbstractMeasureOrDensity

  • parshapes::ValueShapes.AbstractValueShape

  • parbounds::BAT.AbstractVarBounds

Note

Fields parbounds and parbounds do not form part of the stable public API and are subject to change without deprecation.

source
BAT.PriorSubstitutionType
struct PriorSubstitution <: TransformAlgorithm

Substitute 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()
source
BAT.PriorToGaussianType
struct PriorToGaussian <: AbstractTransformTarget

Specifies 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:

  • PriorToGaussian()
source
BAT.PriorToUniformType
struct PriorToUniform <: AbstractTransformTarget

Specifies that posterior densities should be transformed in a way that makes their pior equivalent to a uniform distribution over the unit hypercube.

Constructors:

  • PriorToUniform()
source
BAT.RandResamplingType
struct RandResampling <: AbstractSamplingAlgorithm

Resamples from a given set of samples.

Constructors:

  • RandResampling(; fields...)

Fields:

  • nsamples::Int64: Default: 10 ^ 5
source
BAT.RepetitionWeightingType
struct RepetitionWeighting{T<:AbstractFloat} <: AbstractMCMCWeightingScheme{T}

Sample weighting scheme suitable for sampling algorithms which may repeated samples multiple times in direct succession (e.g. MetropolisHastings). 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()
source
BAT.SampleMedianEstimatorType
struct SampleMedianEstimator <: AbstractMedianEstimator

Get median values from samples using standard Julia statistics functions.

Constructors:

  • SampleMedianEstimator()
source
BAT.SuaveIntegrationType
struct SuaveIntegration <: IntegrationAlgorithm

SuaveIntegration integration algorithm.

Constructors:

  • SuaveIntegration(; fields...)

Fields:

  • trafo::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

Note

This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).

source
BAT.VEGASIntegrationType
struct VEGASIntegration <: IntegrationAlgorithm

VEGASIntegration integration algorithm.

Constructors:

  • VEGASIntegration(; fields...)

Fields:

  • trafo::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

Note

This functionality is only available when the Cuba package is loaded (e.g. via import CUBA).

source
BAT.FixedNBinsType
FixedNBins(nbins::Int)

Selects a fixed number of bins.

Constructor: FixedNBins(; fields...)

Fields:

  • nbins::Int64: Default: 200
source
BAT.ToRealVectorType
struct ToRealVector <: AbstractTransformTarget

Specifies that the input should be transformed into a measure over the space of real-valued flat vectors.

Constructors:

  • ToRealVector()
source
BAT.ConvergenceTestType
abstract type ConvergenceTest

Abstract type for integrated autocorrelation length estimation algorithms.

source