API
Types
BAT.AbstractDensity
BAT.AbstractMCMCCallback
BAT.AbstractPosteriorDensity
BAT.AbstractProposalDist
BAT.AbstractSamplingAlgorithm
BAT.AdaptiveMetropolisTuning
BAT.DataSet
BAT.DistLikeDensity
BAT.GelmanRubinConvergence
BAT.GenericDensity
BAT.HMIData
BAT.HMISettings
BAT.IntegrationVolume
BAT.IntegrationVolume
BAT.MCMCCallbackWrapper
BAT.NamedTupleDist
BAT.OnlineMvMean
BAT.OnlineUvMean
BAT.PointCloud
BAT.PointCloud
BAT.PosteriorSampleVector
BAT.RandSampling
BAT.RandomResampling
BAT.SearchResult
BAT.WhiteningResult
Functions and macros
BAT.apply_bounds
BAT.apply_bounds
BAT.apply_bounds!
BAT.autocrl
BAT.bat_integrate
BAT.bat_read
BAT.bat_read
BAT.bat_sample
BAT.bat_sample
BAT.bat_sampler
BAT.bat_stats
BAT.bat_write
BAT.bat_write
BAT.calculate_localmode
BAT.create_hypercube
BAT.create_hyperrectangle
BAT.default_sampling_algorithm
BAT.density_logval
BAT.distribution_logpdf
BAT.distribution_logpdf!
BAT.drop_low_weight_samples
BAT.effective_sample_size
BAT.effective_sample_size
BAT.effective_sample_size
BAT.eval_density_logval
BAT.eval_prior_posterior_logval!
BAT.eval_prior_posterior_logval_strict!
BAT.find_hypercube_centers
BAT.fromuhc
BAT.fromuhc!
BAT.fromui
BAT.gr_Rsqr
BAT.hm_init
BAT.hm_integrate!
BAT.hm_whiteningtransformation!
BAT.hyperrectangle_creationproccess!
BAT.initial_params!
BAT.issymmetric_around_origin
BAT.log_volume
BAT.modify_hypercube!
BAT.modify_integrationvolume!
BAT.nparams
BAT.nparams
BAT.nparams
BAT.param_bounds
BAT.param_bounds
BAT.params_shape
BAT.proposal_rand!
BAT.spatialvolume
BAT.wgt_effective_sample_size
Documentation
BAT.AbstractDensity
— Type.AbstractDensity
Subtypes of AbstractDensity
only have to imlement the function
BAT.density_logval
However, densities with known parameters bounds should also implement
BAT.param_bounds
If the parameter bounds are unkown, but the number of parameters is known, the function
BAT.nparams
should be implemented directly (usually it is inferred from the bounds).
Densities that support named parameters should also implement
BAT.params_shape
BAT.AbstractMCMCCallback
— Type.AbstractMCMCCallback <: Function
Subtypes (e.g. X
) must support
(::X)(level::Integer, chain::MCMCIterator) => nothing
(::X)(level::Integer, tuner::AbstractMCMCTuner) => nothing
to be compabtible with mcmc_iterate!
, mcmc_tune_burnin!
, etc.
BAT.AbstractPosteriorDensity
— Type.abstract type AbstractPosteriorDensity <: AbstractDensity end
Abstract super-type for posterior probability densities.
BAT.AbstractProposalDist
— Type.AbstractProposalDist
The following functions must be implemented for subtypes:
BAT.distribution_logpdf
BAT.proposal_rand!
BAT.nparams
, returning the number of parameters (i.e. dimensionality).LinearAlgebra.issymmetric
, indicating whether p(a -> b) == p(b -> a) holds true.
In some cases, it may be desirable to override the default implementation of BAT.distribution_logpdf!
.
BAT.AdaptiveMetropolisTuning
— Type.AdaptiveMetropolisTuning(...)
Ajusts the proposal function based on the acceptance ratio and covariance of the previous samples.
BAT.DistLikeDensity
— Type.DistLikeDensity <: AbstractDensity
A density that implements part of the Distributions.Distribution
interface. Such densities are suitable to be used as a priors.
Subtypes of DistLikeDensity
are required to support more functionality than a AbstractDensity
, but less than a Distribution{Multivariate,Continuous}
.
The following functions must be implemented for subtypes:
BAT.density_logval
BAT.param_bounds
BAT.params_shape
Distributions.sampler
Statistics.cov
Prior densities that support named parameters should also implement
BAT.params_shape
A d::Distribution{Multivariate,Continuous}
can be converted into (wrapped in) an DistLikeDensity
via conv(DistLikeDensity, d)
.
BAT.GelmanRubinConvergence
— Type.GelmanRubinConvergence
Gelman-Rubin $$maximum(R^2)$$ convergence test.
BAT.GenericDensity
— Type.GenericDensity{F<:Function} <: AbstractDensity
Constructors:
GenericDensity(log_f)
Turns the logarithmic density function log_f
into a BAT-compatible AbstractDensity
. log_f
must support
`log_f(params::Any)::Real`
It must be safe to execute log_f
in parallel on multiple threads and processes.
BAT.HMIData
— Type.HMIData{T<:AbstractFloat, I<:Integer}
Note: AHMI-internal, not part of stable API.
Includes all the informations of the integration process, including a list of hyper-rectangles, the results of the whitening transformation, the starting ids, and the average number of points and volume of the created hyper-rectangles.
Variables
- 'dataset1' : Data Set 1
- 'dataset2' : Data Set 2
- 'whiteningresult' : contains the whitening matrix and its determinant, required to scale the final integral estimate
- 'volumelist1' : An array of integration volumes created using dataset1, but filled with samples from dataset2
- 'volumelist2' : An array of integration volumes created using dataset2, but filled with samples from dataset1
- 'cubelist1' : An array of small hyper-cubes created around seeding samples of dataset 1
- 'cubelist2' : An array of small hyper-cubes created around seeding samples of dataset 2
- 'iterations1' : The number of volume adapting iterations for the creating volumelist1
- 'iterations2' : The number of volume adapting iterations for the creating volumelist2
- 'rejectedrects1' : An array of ids, indicating which hyper-rectangles of volumelist1 were rejected due to trimming
- 'rejectedrects2' : An array of ids, indicating which hyper-rectangles of volumelist2 were rejected due to trimming
- 'integralestimates' : A dictionary containing the final integral estimates with uncertainty estimation using different uncertainty estimators. Also includes all intermediate results required for the integral estimate combination
BAT.NamedTupleDist
— Type.NamedTupleDist <: MultivariateDistribution
NamedTupleDist <: MultivariateDistribution
A distribution with NamedTuple
-typed variates.
Can be used to describe the distribution of each parameter in a set of named parameters. If the distribution is used as a Bayesian prior, the NamedTupleDist
then specifies the prior on each named parameter.
A NamedTupleDist
implies a NamedTupleShape
:
valshape(d::NamedTupleDist) isa NamedTupleShape
BAT.OnlineMvMean
— Type.OnlineMvMean{T<:AbstractFloat} <: AbstractVector{T}
Multi-variate mean implemented via Kahan-Babuška-Neumaier summation.
BAT.OnlineUvMean
— Type.OnlineUvMean{T<:AbstractFloat}
Univariate mean implemented via Kahan-Babuška-Neumaier summation.
BAT.PosteriorSampleVector
— Type.PosteriorSampleVector
Type alias for StructArrays.StructArray{<:PosteriorSample,...}
.
BAT.bat_integrate
— Function.bat_integrate(
posterior::BAT.AnyPosterior,
)::PosteriorSampleVector
Calculate the integral (evidence) of posterior
.
Returns a NamedTuple: (integral = x::Measurement.Measurement, ...)
Result properties not listed here are algorithm-specific and are not part of the stable BAT API.
posterior
may be a
Distributions.MultivariateDistribution
Uses the AHMI algorithm by default.
BAT.bat_read
— Function.bat_read(filename::AbstractString, data)
bat_read(fn_with_subpath::Tuple{AbstractString, AbstractString}, data)
Read data from a file filename
, resp. from an internal sub-path of the file (if supported by the file format), e.g. an HDF5 group.
Currently supported file formats are:
- HDF5 with file extension ".h5" or ".hdf5"
BAT.bat_read
— Method.bat_read(src::HDF5.DataFile)
bat_read(src_with_subpath::Tuple{HDF5.DataFile, AbstractString})
Read data from HDF5 file or group src
(optionally from an HDF5-path relative to src
).
BAT.bat_sample
— Function.bat_sample(
[rng::AbstractRNG],
posterior::BAT.AnyPosterior,
n::BAT.AnyNSamples,
[algorithm::BAT.AbstractSamplingAlgorithm]
)::PosteriorSampleVector
Draw n
samples from posterior
.
Returns a NamedTuple of the shape
(
samples = X::PosteriorSampleVector,...
stats = s::@test stats isa NamedTuple{(:mode,:mean,:cov,...)},
...
)
Result properties not listed here are algorithm-specific and are not part of the stable BAT API.
posterior
may be a
Distributions.MultivariateDistribution
Depending on the type of posterior, n
may be of type
Integer
: Number of samplesTuple{Integer,Integer}
: Tuple of number of samples per sample source and number of sample sources (e.g. number of MCMC chains). The total number of samples isproduct(n)
.
Depending on the type of posterior
, the number of samples returned may be somewhat larger or smaller than specified by product(n)
.
Also depending on the posterior
type, the samples may be independent or correlated (e.g. when using MCMC).
BAT.bat_sample
— Method.function bat_sample(
rng::AbstractRNG,
posterior::AbstractPosteriorDensity,
n::Union{Integer,Tuple{Integer,Integer}},
algorithm::MCMCAlgorithm;
max_nsteps::Integer,
max_time::Real,
tuning::AbstractMCMCTunerConfig,
init::MCMCInitStrategy,
burnin::MCMCBurninStrategy,
convergence::MCMCConvergenceTest,
strict::Bool = false,
filter::Bool = true
)
Sample posterior
via Markov chain Monte Carlo (MCMC).
n
must be either a tuple (nsteps, nchains)
or an integer. nchains
specifies the (approximate) number of MCMC steps per chain, nchains
the number of MCMC chains. If n is an integer, it is interpreted as nsteps * nchains
, and the number of steps and chains are chosen automatically.
BAT.bat_sampler
— Function.bat_sampler(d::Distribution)
Tries to return a BAT-compatible sampler for Distribution d. A sampler is BAT-compatible if it supports random number generation using an arbitrary AbstractRNG
:
rand(rng::AbstractRNG, s::SamplerType)
rand!(rng::AbstractRNG, s::SamplerType, x::AbstractArray)
If no specific method of bat_sampler
is defined for the type of d
, it will default to sampler(d)
, which may or may not return a BAT-compatible sampler.
BAT.bat_stats
— Method.bat_stats(samples::PosteriorSampleVector)
Calculated parameter statistics on samples
. Returns a NamedTuple{(:mode,:mean,:cov,...)}
. Result properties not listed here are not part of the stable BAT API and subject to change.
BAT.bat_write
— Function.bat_write(filename::AbstractString, data)
bat_write(fn_with_subpath::Tuple{AbstractString, AbstractString}, data)
Write data to a file filename
, resp. to an internal sub-path of the file (if supported by the file format), e.g. an HDF5 group.
Currently supported file formats are:
- HDF5 with file extension ".h5" or ".hdf5"
BAT.bat_write
— Method.bat_write(dest::HDF5.DataFile, data)
bat_write(dest_with_subpath::Tuple{HDF5.DataFile, AbstractString}, data)
Write data
to HDF5 file or group dest
(optionally to an HDF5-path relative to dest
).
data
must be a table (i.e. implement the Tables.jl API).
BAT.density_logval
— Function.density_logval(density::AbstractDensity, params::Any)
Compute log of value of a multi-variate density function at the given parameter values.
Input:
density
: density functionparams
: parameter values
Note: If density_logval
is called with out-of-bounds parameters (see param_bounds
), the behaviour is undefined. The result for parameters that are not within bounds is implicitly -Inf
, but it is the caller's responsibility to handle these cases.
BAT.distribution_logpdf
— Function.distribution_logpdf(
pdist::AbstractProposalDist,
params_new::AbstractVector,
params_old:::AbstractVector
)
Analog to distribution_logpdf!
, but for a single parameter vector.
BAT.distribution_logpdf!
— Function.distribution_logpdf!(
p::AbstractArray,
pdist::AbstractProposalDist,
params_new::Union{AbstractVector,VectorOfSimilarVectors},
params_old:::Union{AbstractVector,VectorOfSimilarVectors}
)
Returns log(PDF) value of pdist
for transitioning from old to new parameter values for multiple parameter sets.
end
Input:
params_new
: New parameter values (column vectors)params_old
: Old parameter values (column vectors)
Output is stored in
p
: Array of PDF values, length must match, shape is ignored
Array size requirements:
size(params_old, 1) == size(params_new, 1) == length(pdist)
size(params_old, 2) == size(params_new, 2)
orsize(params_old, 2) == 1
size(params_new, 2) == length(p)
Implementations of distribution_logpdf!
must be thread-safe.
BAT.fromuhc!
— Function.fromuhc!(Y::AbstractVector, X::AbstractVector, vol::SpatialVolume)
fromuhc!(Y::VectorOfSimilarVectors, X::VectorOfSimilarVectors, vol::SpatialVolume)
Bijective transformation of coordinates X
within the unit hypercube to coordinates Y
in vol
. If X
and Y
are matrices, the transformation is applied to the column vectors. Use Y === X
to transform in-place.
Use inv(fromuhc!)
to get the the inverse transformation.
BAT.fromuhc
— Method.fromuhc(X::AbstractVector, vol::SpatialVolume)
fromuhc(X::VectorOfSimilarVectors, vol::SpatialVolume)
Bijective transformation from unit hypercube to vol
. See fromuhc!
.
Use inv(fromuhc)
to get the the inverse transformation.
BAT.fromui
— Function.y = fromui(x::Real, lo::Real, hi::Real)
y = fromui(x::Real, lo_hi::ClosedInterval{<:Real})
Linear bijective transformation from the unit inverval (i.e. x ∈ 0..1
) to y ∈ lo..hi
.
Use inv(fromui)
to get the the inverse transformation.
Use @inbounds
to disable range checking on the input value.
BAT.issymmetric_around_origin
— Function.issymmetric_around_origin(d::Distribution)
Returns true
(resp. false
) if the Distribution is symmetric (resp. non-symmetric) around the origin.
BAT.log_volume
— Function.log_volume(vol::SpatialVolume)
Get the logarithm of the volume of the space in vol
.
BAT.nparams
— Function.nparams(X::Union{AbstractParamBounds,MCMCIterator,...})
Get the number of parameters of X
.
BAT.nparams
— Method.nparams(density::AbstractDensity)::Union{Int,Missing}
Get the number of parameters of density
. May return missing
, if the density supports a variable number of parameters.
BAT.nparams
— Method.nparams(density::DistLikeDensity)::Int
Get the number of parameters of prior density density
. Must not be missing
, prior densities must have a fixed number of parameters. By default, the number of parameters is inferred from the parameter bounds.
BAT.param_bounds
— Function.param_bounds(density::DistLikeDensity)::AbstractParamBounds
Get the parameter bounds of density
. Must not be missing
.
BAT.param_bounds
— Method.param_bounds(
density::AbstractDensity
)::Union{AbstractParamBounds,Missing}
Get the parameter bounds of density
. See density_logval
for the implications and handling of the bounds. If the bounds are missing, density_logval
must be prepared to handle any parameter values.
BAT.params_shape
— Function.params_shape(
density::AbstractDensity
)::Union{ValueShapes.AbstractValueShape,Missing}
params_shape(
density::DistLikeDensity
)::ValueShapes.AbstractValueShape
params_shape(
distribution::Distributions.Distribution
)::ValueShapes.AbstractValueShape
Get the shapes of parameters of density
.
For prior densities, the result must not be missing
, but may be nothing
if the prior only supports flat parameter vectors.
BAT.proposal_rand!
— Function.function proposal_rand!(
rng::AbstractRNG,
pdist::GenericProposalDist,
params_new::Union{AbstractVector,VectorOfSimilarVectors},
params_old::Union{AbstractVector,VectorOfSimilarVectors}
)
Generate one or multiple proposed parameter vectors, based on one or multiple previous parameter vectors.
Input:
rng
: Random number generator to usepdist
: Proposal distribution to useparams_old
: Old parameter values (vector or column vectors, if a matrix)
Output is stored in
params_new
: New parameter values (vector or column vectors, if a matrix)
The caller must guarantee:
size(params_old, 1) == size(params_new, 1)
size(params_old, 2) == size(params_new, 2)
orsize(params_old, 2) == 1
params_new !== params_old
(no aliasing)
Implementations of proposal_rand!
must be thread-safe.
BAT.spatialvolume
— Function.spatialvolume(b::ParamVolumeBounds)::SpatialVolume
Returns the spatial volume that defines the parameter bounds.
BAT.AbstractSamplingAlgorithm
— Type.BAT.AbstractSamplingAlgorithm
Abstract type for BAT sampling algorithms. See bat_sample
.
BAT.DataSet
— Type.DataSet{T<:AbstractFloat, I<:Integer}
Note: AHMI-internal, not part of stable API.
Holds the MCMC output. For construction use constructor: function DataSet{T<:Real}(data::Matrix{T}, logprob::Vector{T}, weights::Vector{T})
Variables
- 'data' : An P x N array with N data points with P parameters.
- 'logprob' : The logarithmic probability for each samples stored in an array
- 'weights' : How often each sample occurred. Set to an array of ones if working directly on MCMC output
- 'ids' : Array which is used to assign each sample to a batch, required for the cov. weighed uncertainty estimation
- .sortids : an array of indices which stores the original ordering of the samples (the space partitioning tree reorders the samples), required to calculate an effective sample size.
- 'N' : number of samples
- 'P' : number of parameters
- 'nsubsets' : the number of batches
- 'iswhitened' : a boolean value which indicates whether the data set is iswhitened
- 'isnew' : a boolean value which indicates whether the data set was swapped out with a different one (it is possible to redo the integration with a different sample set using previously generated hyper-rectangles)
- 'partitioningtree' : The space partitioning tree, used to efficiently identify samples in a point cloud
- 'startingIDs' : The Hyper-Rectangle Seed Samples are stored in this array
- 'tolerance' : A threshold required for the hyper-rectangle creation process.
BAT.HMISettings
— Type.HMISettings
Note: AHMI-internal, not part of stable API.
holds the settings for the hm_integrate function. There are several default constructors available: HMIFastSettings() HMIStandardSettings() HMIPrecisionSettings()
#Variables
- 'whitening_method::Symbol' : which whitening method to use
- 'max_startingIDs::Integer' : influences how many starting ids are allowed to be generated
- 'maxstartingIDsfraction::AbstractFloat' : how many points are considered as possible starting points as a fraction of total points available
- 'rect_increase::AbstractFloat' : describes the procentual rectangle volume increase/decrease during hyperrectangle creation. Low values can increase the precision if enough points are available but can cause systematically wrong results if not enough points are available.
- 'useallrects::Bool' : All rectangles are used for the integration process no matter how big their overlap is. If enabled the rectangles are weighted by their overlap.
- 'useMultiThreading' : activate multithreading support.
- 'warning_minstartingids' : the required minimum amount of starting samples
- 'dotrimming' : determines whether the integral estimates are trimmed (1σ trim) before combining them into a final result (more robust)
- 'uncertaintyestimators' : A dictionary of different uncertainty estimator functions. Currently three functions are available: hmcombineresultslegacy! (outdated, overestimates uncertainty significantly in higher dimensions), hmcombineresultscovweighted! (very fast) and hmcombineresults_analyticestimation! (recommended)
end
BAT.IntegrationVolume
— Type.IntegrationVolume{T<:AbstractFloat, I<:Integer}
Note: AHMI-internal, not part of stable API.
Variables
- 'pointcloud' : holds the point cloud of the integration volume
- 'spatialvolume' : the boundaries of the integration volume
- 'volume' : the volume
Hold the point cloud and the spatial volume for integration.
BAT.IntegrationVolume
— Method.IntegrationVolume(dataset::DataSet{T, I}, spvol::HyperRectVolume{T}, searchpts::Bool = true)::IntegrationVolume{T, I}
Note: AHMI-internal, not part of stable API.
creates an integration region by calculating the point cloud an the volume of the spatial volume.
BAT.MCMCCallbackWrapper
— Type.MCMCCallbackWrapper{F} <: AbstractMCMCCallback
Wraps a callable object to turn it into an AbstractMCMCCallback
.
Constructor:
MCMCCallbackWrapper(f::Any)
f
needs to support the call syntax of an AbstractMCMCCallback
.
BAT.PointCloud
— Type.PointCloud{T<:AbstractFloat, I<:Integer}
Note: AHMI-internal, not part of stable API.
Stores the information of the points of an e.g. HyperRectVolume
Variables
- 'maxLogProb' : The maximum log. probability of one of the points inside the hyper-rectangle
- 'minLogProb' : The minimum log. probability of one of the points inside the hyper-rectangle
- 'maxWeightProb' : the weighted max. log. probability
- 'minWeightProb' : the weighted min. log. probability
- 'probfactor' : The probability factor of the hyper-rectangle
- 'probweightfactor' : The weighted probability factor
- 'points' : The number of points inside the hyper-rectangle
- 'pointIDs' : the IDs of the points inside the hyper-rectangle, might be empty because it is optional and costs performance
- 'searchres' : used to boost performance
BAT.PointCloud
— Method.PointCloud{T<:AbstractFloat, I<:Integer}(dataset::DataSet{T, I}, hyperrect::HyperRectVolume{T}, searchpts::Bool = false)::PointCloud
Note: AHMI-internal, not part of stable API.
creates a point cloud by searching the data tree for points which are inside the hyper-rectangle The parameter searchpts determines if an array of the point IDs is created as well
BAT.RandSampling
— Type.BAT.RandSampling
Constructors:
BAT.RandSampling()
Sample via Random.rand
. Only supported for posteriors of type Distributions.MultivariateDistribution
and BAT.DistLikeDensity
.
BAT.RandomResampling
— Type.BAT.RandomResampling
Constructors:
BAT.RandomResampling()
Resample from a given set of samples.
BAT.SearchResult
— Type.SearchResult{T<:AbstractFloat, I<:Integer}
Note: AHMI-internal, not part of stable API.
Stores the results of the space partitioning tree's search function
Variables
- 'pointIDs' : the IDs of samples found, might be empty because it is optional
- 'points' : The number of points found.
- 'maxLogProb' : the maximum log. probability of the points found.
- 'minLogProb' : the minimum log. probability of the points found.
- 'maxWeightProb' : the weighted minimum log. probability found.
- 'minWeightProb' : the weighted maximum log. probfactor found.
BAT.WhiteningResult
— Type.WhiteningResult{T<:AbstractFloat}
Note: AHMI-internal, not part of stable API.
Stores the information obtained during the Whitening Process
Variables
- 'determinant' : The determinant of the whitening matrix
- 'targetprobfactor' : The suggested target probability factor
- 'whiteningmatrix' : The whitening matrix
- 'meanvalue' : the mean vector of the input data
BAT.apply_bounds
— Function.apply_bounds(x::Real, interval::ClosedInterval, boundary_type::BoundsType)
Specify lower and upper bound via interval
.
BAT.apply_bounds!
— Function.apply_bounds!(params::AbstractVector, bounds::AbstractParamBounds)
Apply bounds
to parameters params
.
BAT.apply_bounds
— Method.apply_bounds(x::<:Real, lo::<:Real, hi::<:Real, boundary_type::BoundsType)
Apply lower/upper bound lo
/hi
to value x
. boundary_type
may be hard_bounds
, cyclic_bounds
or reflective_bounds
.
BAT.autocrl
— Method.autocrl(xv::AbstractVector{T}, kv::AbstractVector{Int} = Vector{Int}())
autocorrelation := Σ Cov[xi,x(i+k)]/Var[x]
Computes the autocorrelations at various leg k of the input vector (time series) xv. The vector kv is the collections of lags to take into account
BAT.calculate_localmode
— Method.calculate_localmode(hist)
Calculates the modes of a 1d statsbase histogram. A vector containing the bin-center(s) of the heighest bin(s) is returned.
BAT.create_hypercube
— Method.create_hypercube{T<:Real}(origin::Vector{T}, edgelength::T)::HyperRectVolume
Note: AHMI-internal, not part of stable API.
creates a hypercube shaped spatial volume
BAT.create_hyperrectangle
— Method.create_hyperrectangle(...)
Note: AHMI-internal, not part of stable API.
This function creates a hyper-rectangle around each starting sample. It starts by building a hyper-cube and subsequently adapts each face individually, thus turning the hyper-cube into a hyper-rectangle. The faces are adjusted in a way to match the shape of the distribution as best as possible.
BAT.default_sampling_algorithm
— Function.BAT.default_sampling_algorithm(posterior)
Get BAT's default sampling algorithm for posterior
.
BAT.drop_low_weight_samples
— Function.drop_low_weight_samples(
samples::PosteriorSampleVector,
fraction::Real = 10^-4
)
Drop fraction
of the total probability mass from samples to filter out the samples with the lowest weight.
Note: BAT-internal function, not part of stable API.
BAT.effective_sample_size
— Method.effective_sample_size(params::AbstractArray, weights::AbstractVector; with_weights=true)
Effective size estimation for a (multidimensional) ElasticArray. By default applies the Kish approximation with the weigths available, but can be turned off (with_weights=false).
BAT.effective_sample_size
— Method.effective_sample_size(samples::PosteriorSampleVector; with_weights=true)
Effective size estimation for a (multidimensional) PosteriorSampleVector. By default applies the Kish approximation with the weigths available, but can be turned off (with_weights=false).
BAT.effective_sample_size
— Method.Effective size estimation for a vector of samples xv. If a weight vector w is provided, the Kish approximation is applied.
By default computes the autocorrelation up to the square root of the number of entries in the vector, unless an explicit list of lags is provided (kv).
BAT.eval_density_logval
— Method.eval_density_logval(
density::AbstractDensity,
params::AbstractVector{<:Real},
parshapes::ValueShapes.AbstractValueShape
)
Internal function to evaluate density log-value, calls density_logval
.
parshapes
must be compatible with params_shape(density)
.
Checks that:
- The number of parameters of
density
(if known) matches the length ofparams
. - The return value of
density_logval
is notNaN
. - The return value of
density_logval
is less than+Inf
.
BAT.eval_prior_posterior_logval!
— Function.BAT.eval_prior_posterior_logval!(
T::Type{<:Real},
density::AbstractDensity,
params::AbstractVector{<:Real}
)
First apply bounds to the parameters, compute prior and posterior log values by via eval_density_logval
.
May modify params
.
Returns a NamedTuple{(:log_prior, :log_posterior),Tuple{T,T}}
Guarantees that :
- If parameters are still out of bounds after applying bounds,
density_logval
is not called for either prior or likelihood. - If
density_logval
for prior returns-Inf
,density_logval
is not called for likelihood.
In both cases, T(-Inf)
is returned for both prior and posterior.
BAT.eval_prior_posterior_logval_strict!
— Method.BAT.eval_prior_posterior_logval_strict!(
density::AbstractDensity,
params::AbstractVector{<:Real}
)
First apply bounds to the parameters, compute prior and posterior log values by via eval_density_logval
.
May modify params
.
Returns a NamedTuple{(:log_prior, :log_posterior),Tuple{T,T}}
. T is inferred from value returned by eval_density_logval
for the likelihood.
Guarantees that :
- If parameters are still out of bounds after applying bounds,
density_logval
is not called for either prior or likelihood. - If
density_logval
for prior returns-Inf
,density_logval
is not called for likelihood.
In both cases, an exception is thrown.
BAT.find_hypercube_centers
— Method.find_hypercube_centers(dataset::DataSet{T, I}, whiteningresult::WhiteningResult, settings::HMISettings)::Vector{I}
Note: AHMI-internal, not part of stable API.
finds possible starting points for the hyperrectangle creation
BAT.gr_Rsqr
— Method.gr_Rsqr(stats::AbstractVector{<:MCMCBasicStats})
Gelman-Rubin $$R^2$$ for all parameters.
BAT.hm_init
— Method.hm_init!(result, settings)
Note: AHMI-internal, not part of stable API.
Sets the global multithreading setting and ensures that a minimum number of samples, dependent on the number of dimensions, are provided.
BAT.hm_integrate!
— Method.hm_integrate!(result, settings = HMIPrecisionSettings())
Note: AHMI-internal, not part of stable API.
This function starts the adaptive harmonic mean integration. See arXiv:1808.08051 for more details. It needs a HMIData struct as input, which holds the samples, in form of a dataset, the integration volumes and other properties, required for the integration, and the final result.
BAT.hm_whiteningtransformation!
— Method.hm_whiteningtransformation!(result, settings)
Note: AHMI-internal, not part of stable API.
Applies a whitening transformation to the samples. A custom whitening method can be used by overriding settings.whitening_function!
BAT.hyperrectangle_creationproccess!
— Method.hyperrectangle_creationproccess!(...)
Note: AHMI-internal, not part of stable API.
This function assigns each thread its own hyper-rectangle to build, if in multithreading-mode.
BAT.initial_params!
— Function.BAT.initial_params!(
params::Union{AbstractVector{<:Real},VectorOfSimilarVectors{<:Real}},
rng::AbstractRNG,
posterior::AbstractPosteriorDensity,
algorithm::MCMCAlgorithm
)::typeof(params)
Fill params
with random initial parameters suitable for posterior
and algorithm
. The default implementation will try to draw the initial parameters from the prior of the posterior.
BAT.modify_hypercube!
— Method.create_hypercube!{T<:Real}(origin::Vector{T}, edgelength::T)::HyperRectVolume
Note: AHMI-internal, not part of stable API.
resizes a hypercube shaped spatial volume
BAT.modify_integrationvolume!
— Method.modify_integrationvolume!(intvol::IntegrationVolume{T, I}, dataset::DataSet{T, I}, spvol::HyperRectVolume{T}, searchpts::Bool = true)
Note: AHMI-internal, not part of stable API.
updates an integration volume with new boundaries. Recalculates the pointcloud and volume.
BAT.wgt_effective_sample_size
— Method.wgt_effective_sample_size(w::AbstractVector{T})
Kish's approximation for weighted samples effectivesamplesize estimation. Computes the weighting factor for weigthed samples, where w is the vector of weigths.