Data Types
Problem & Model
The BolfiProblem
structure contains all information about the inference problem, as well as the model hyperparameters.
BOLFI.BolfiProblem
— TypeBolfiProblem(X, Y; kwargs...)
BolfiProblem(data::ExperimentData; kwargs...)
BolfiProblem(::BossProblem, ::Likelihood, x_prior::MultivariateDistribution)
BolfiProblem(::BossProblem, ::Likelihood, x_prior::MultivariateDistribution, y_sets::Matrix{Bool})
Defines the likelihood-free inference problem and stores all data.
Args
The initial data are provided either as two column-wise matrices X
and Y
with inputs and outputs of the simulator respectively, or as an instance of BOSS.ExperimentData
.
Currently, at least one datapoint has to be provided (purely for implementation reasons).
Kwargs
f::Any
: The simulation to be queried for data.bounds::AbstractBounds
: The basic box-constraints onx
. This field is mandatory.discrete::AbstractVector{<:Bool}
: Can be used to designate some dimensions of the domain as discrete.cons::Union{Nothing, Function}
: Used to define arbitrary nonlinear constraints onx
. Feasible pointsx
must satisfyall(cons(x) .> 0.)
. An appropriate acquisition maximizer which can handle nonlinear constraints must be used ifcons
is provided. (SeeBOSS.AcquisitionMaximizer
.)acquisition::BolfiAcquisition
: Defines the acquisition function.kernel::Kernel
: The kernel used in the GP. Defaults to theMatern32Kernel()
.lengthscale_priors::AbstractVector{<:MultivariateDistribution}
: The prior distributions for the length scales of the GP. Thelengthscale_priors
should be a vector ofy_dim
x_dim
-variate distributions wherex_dim
andy_dim
are the dimensions of the input and output of the model respectively.amplitude_priors::AbstractVector{<:UnivariateDistribution}
: The prior distributions for the amplitude hyperparameters of the GP. Theamplitude_priors
should be a vector ofy_dim
univariate distributions.noise_std_priors::AbstractVector{<:UnivariateDistribution}
: The prior distributions of the standard deviations the Gaussian simulator noise on each dimension of the outputy
.likelihood::Likelihood
: The likelihood of the experiment observationy_o
.x_prior::MultivariateDistribution
: The priorp(x)
on the input parameters.y_sets::Matrix{Bool}
: Optional parameter intended for advanced usage. The binary columns define subsetsy_1, ..., y_m
of the observation dimensions withiny
. The algorithm then trains multiple posteriorsp(θ|y_1), ..., p(θ|y_m)
simultaneously. The posteriors can be compared after the run is completed to see which observation subsets are most informative.
Likelihood
The abstract type Likelihood
represents the likelihood distribution of the observation y_o
.
BOLFI.Likelihood
— TypeRepresents the assumed likelihood of the experiment observation $y_o$.
Defining a Custom Likelihood
To define a custom likelihood, create a new struct CustomLike <: Likelihood
and implement the following methods;
loglike(::CustomLike, z::AbstractVector{<:Real})
wherez
is the simulator outputapprox_likelihood(::CustomLike, bolfi, gp_post)
: Necessary to be able to useapprox_posterior
.likelihood_mean(::CustomLike, bolfi, gp_post)
: Necessary to be able to useposterior_mean
.sq_likelihood_mean(::CustomLike, bolfi, gp_post)
: Necessary to be able to useposterior_variance
.get_subset(::CustomLike, y_set::AsbtractVector{<:Bool})
: Necessary forBolfiProblem
s where!isnothing(problem.y_sets)
.
The NormalLikelihood
assumes that the observation y_o
has been drawn from a Gaussian distribution with a known diagonal covariance matrix with the std_obs
values on the diagonal. The simulator is used to learn the mean function.
BOLFI.NormalLikelihood
— TypeNormalLikelihood(; y_obs, std_obs)
The observation is assumed to have been generated from a normal distribution as y_o \sim Normal(f(x), Diagonal(std_obs))
. We can use the simulator to query z = f(x)
.
Kwargs
y_obs::Vector{Float64}
: The observed values from the real experiment.std_obs::Union{Vector{Float64}, Nothing}
: The standard deviations of the Gaussian observation noise on each dimension of the "ground truth" observation. (If the observation is considered to be generated from the simulator and not some "real" experiment, providestd_obs = nothing
` and the adaptively trained simulation noise deviation will be used in place of the experiment noise deviation as well. This may be the case for some toy problems or benchmarks.)
The LogNormalLikelihood
assumes that the observation y_o
has been drawn from a log-normal distribution with a known diagonal covariance matrix with the std_obs
values on the diagonal. The simulator is used to learn the mean function.
BOLFI.LogNormalLikelihood
— TypeLogNormalLikelihood(; y_obs, std_obs)
The observation is assumed to have been generated from a normal distribution as y_o \sim LogNormal(f(x), Diagonal(std_obs))
. We can use the simulator to query z = f(x)
.
In many cases, one may want to take the logarithm of the output of the simulator. Meaning, if one has simulator z = sim(x)
, one would define f
as y = f(x) = log(sim(x))
. This way, the y
values with high likelihood will have similar values to the z
values.
Kwargs
y_obs::Vector{Float64}
: The observed values from the real experiment.std_obs::Vector{Float64}
: The standard deviations of the LogNormal observation noise. (If the observation is considered to be generated from the simulator and not some "real" experiment, providestd_obs = nothing
` and the adaptively trained simulation noise deviation will be used in place of the experiment noise deviation as well. This may be the case for some toy problems or benchmarks.)
The BinomialLikelihood
assumes that the observation y_o
has been drawn from a Binomial distribution with a known number trials
. The simulator is used to learn the probability parameter p
as a function of the input parameters. The expectation over this likelihood (in case one wants to use posterior_mean
and/or posterior_variance
) is calculated via simple numerical integration on a predefined grid.
BOLFI.BinomialLikelihood
— TypeBinomialLikelihood(; y_obs, trials, kwargs...)
The observation is assumed to have been generated from a Binomial distribution as y_o \sim Binomial(trials, f(x))
. We can use the simulator to query z = f(x)
.
The simulator should only return values between 0 and 1. The GP estimates are clamped to this range.
Kwargs
y_obs::Vector{Int64}
: The observed values from the real experiment.trials::Vector{Int64}
: The number of trials for each observation dimension.int_grid_size::Int64
: The number of samples used to approximate the expected likelihood.
The GutmannNormalLikelihood
is implemented according to the equations of Gutmann et al. in [1,2]. It is defined in a different way to the other Likelihood
s. It is provided mainly for a comparison with NormalLikelihood
.
BOLFI.GutmannNormalLikelihood
— TypeGutmannNormalLikelihood(; std_obs::Vector{Float64})
The observation is assumed to have been generated from a normal distribution as y_o \sim Normal(f(x), Diagonal(std_obs))
. We can use the simulator to query z = f(x)
.
This likelihood follows the equations from Gutmann et al. [1,2].
It is defined in a slightly different way than the other Likelihood
s. Namely, this likelihood requires the simulator to return a single non-negative scalar value, describing the discrepancy \delta
of the simulator outcome from the real observation.
(See that the real observation y_obs
is not provided as a parameter to GutmannNormalLikelihood
as in the case of other Likelihood
s. Instead, the observation y_obs
should be used to calculate the discrepancy \delta
in the simulator f
provided to BolfiProblem
. The simulator should return the discrepancy as a vector [\delta]
of length 1.)
The likelihood is then defined as $P[\delta < \epslion]$ instead of $P[\delta = 0]$, which would correspond more closely to the other Likelihood
s. The parameter \epsilon
is the acceptance threshold.
Kwargs
ϵ::Float64
: The threshold for the discrepancy from the real observation.
References
[1] Gutmann, Michael U., and Jukka Cor. "Bayesian optimization for likelihood-free inference of simulator-based statistical models." Journal of Machine Learning Research 17.125 (2016): 1-47.
[2] Järvenpää, Marko, et al. "Efficient acquisition rules for model-based approximate Bayesian computation." (2019): 595-622.
Acquisition Function
The abstract type BolfiAcquisition
represents the acquisition function.
PostVarAcq, MWMVAcq, InfoGain
BOLFI.BolfiAcquisition
— TypeAn abstract type for BOLFI acquisition functions.
Required API for subtypes of BolfiAcquisition
:
- Implement method
(::CustomAcq)(::Type{<:UniFittedParams}, ::BolfiProblem, ::BolfiOptions) -> (x -> ::Real)
.
Optional API for subtypes of BolfiAcquisition
:
- Implement method
(::CustomAcq)(::Type{<:MultiFittedParams}, ::BolfiProblem, ::BolfiOptions) -> (x -> ::Real)
. A default fallback is provided forMultiFittedParams
, which averages individual acquisition functions for each sample.
The PostVarAcq
can be used to solve LFI problems. It maximizes the posterior variance to select the next evaluation point.
BOLFI.PostVarAcq
— TypePostVarAcq()
Selects the new evaluation point by maximizing the variance of the posterior approximation.
The MWMVAcq
can be used to solve LFSS problems. It maximizes the "mass-weighted mean variance" of the posteriors given by the different sensor sets.
BOLFI.MWMVAcq
— TypeMWMVAcq(; kwargs...)
The Mass-Weighted Mean Variance acquisition function.
Selects the next evaluation point by maximizing a weighted average of the variances of the individual posterior approximations given by different sensor sets. The weights are determined as the total probability mass of the current data w.r.t. each approximate posterior.
Keywords
samples::Int
: The number of samples used to estimate the evidence.
Termination Condition
The abstract type BolfiTermCond
represents the termination condition for the whole BOLFI procedure. Additionally, any BOSS.TermCond
from the BOSS.jl package can be used with BOLFI.jl as well, and it will be automatically converted to a BolfiTermCond
.
BOLFI.BolfiTermCond
— TypeAn abstract type for BOLFI termination conditions.
Implementing custom termination condition:
- Create struct
CustomTermCond <: BolfiTermCond
- Implement method
(::CustomTermCond)(::BolfiProblem) -> ::Bool
The most basic termination condition is the BOSS.IterLimit
, which can be used to simply terminate the procedure after a predefined number of iterations.
BOLFI.jl provides two specialized termination conditions; the AEConfidence
, and the UBLBConfidence
. Both of them estimate the degree of convergence by comparing confidence regions given by two different approximations of the posterior.
BOLFI.AEConfidence
— TypeAEConfidence(; kwargs...)
Calculates the q
-confidence region of the expected and the approximate posteriors. Terminates after the IoU of the two confidence regions surpasses r
.
Keywords
max_iters::Union{Nothing, <:Int}
: The maximum number of iterations.samples::Int
: The number of samples used to approximate the confidence regions and their IoU ratio. Only has an effect ifisnothing(xs)
.xs::Union{Nothing, <:AbstractMatrix{<:Real}}
: Can be used to provide a pre-sampled set of parameter samples from thex_prior
defined inBolfiProblem
.q::Float64
: The confidence value of the confidence regions. Defaults toq = 0.95
.r::Float64
: The algorithm terminates once the IoU ratio surpassesr
. Defaults tor = 0.95
.
BOLFI.UBLBConfidence
— TypeUBLBConfidence(; kwargs...)
Calculates the q
-confidence region of the UB and LB approximate posterior. Terminates after the IoU of the two confidence intervals surpasses r
. The UB and LB confidence intervals are calculated using the GP mean +- n
GP stds.
Keywords
max_iters::Union{Nothing, <:Int}
: The maximum number of iterations.samples::Int
: The number of samples used to approximate the confidence regions and their IoU ratio. Only has an effect ifisnothing(xs)
.xs::Union{Nothing, <:AbstractMatrix{<:Real}}
: Can be used to provide a pre-sampled set of parameter samples from thex_prior
defined inBolfiProblem
.n::Float64
: The number of predictive deviations added/substracted from the GP mean to get the two posterior approximations. Defaults ton = 1.
.q::Float64
: The confidence value of the confidence regions. Defaults toq = 0.8
.r::Float64
: The algorithm terminates once the IoU ratio surpassesr
. Defaults tor = 0.8
.
Miscellaneous
The BolfiOptions
structure can be used to define miscellaneous settings of BOLFI.jl.
BOLFI.BolfiOptions
— TypeBolfiOptions(; kwargs...)
Stores miscellaneous settings.
Keywords
info::Bool
: Settinginfo=false
silences the algorithm.debug::Bool
: Setdebug=true
to print stactraces of caught optimization errors.parallel_evals::Symbol
: Possible values::serial
,:parallel
,:distributed
. Defaults to:parallel
. Determines whether to run multiple objective function evaluations within one batch in serial, parallel, or distributed fashion. (Only has an effect if batching AM is used.)callback::Union{<:BossCallback, <:BolfiCallback}
: If provided, the callback will be called before the BO procedure starts and after every iteration.
The abstract type BolfiCallback
can be derived to define a custom callback, which will be called once before the BOLFI procedure starts, and subsequently in every iteration.
For an example usage of this functionality, see the example in the package repository, where a custom callback is used to create the plots.
BOLFI.BolfiCallback
— TypeIf a callback cb
of type BolfiCallback
is defined in BolfiOptions
, the method cb(::BolfiProblem; kwargs...)
will be called in every iteration.
cb(problem::BolfiProblem;
model_fitter::BOSS.ModelFitter,
acq_maximizer::BOSS.AcquisitionMaximizer,
acquisition::BOLFI.AcqWrapper, # `BolfiAcquisition` wrapped into `AcqWrapper`
term_cond::TermCond, # either `BOSS.TermCond` or a `BolfiTermCond` wrapped into `TermCondWrapper`
options::BossOptions,
first::Bool,
)
References
[1] Gutmann, Michael U., and Jukka Cor. "Bayesian optimization for likelihood-free inference of simulator-based statistical models." Journal of Machine Learning Research 17.125 (2016): 1-47.
[2] Järvenpää, Marko, et al. "Efficient acquisition rules for model-based approximate Bayesian computation." (2019): 595-622.