API

SimulatedAnnealingABC.sabcFunction
sabc(f_dist::Function, prior::Distribution, args...;
      n_particles = 100, n_simulation = 10_000,
      algorithm = :single_eps,
      propsal =  DifferentialEvolution(γ0=2.38/sqrt(2*n_para)).
      resample = 2*n_particles,
      v=1.0, δ=0.1,
      checkpoint_history = 1,
      show_progressbar::Bool = !is_logging(stderr),
      show_checkpoint = is_logging(stderr) ? 100 : Inf,
      kwargs...)

Simulated Annealing Approximate Bayesian Inference Algorithm

Arguments

  • f_dist: Function that returns one or more distances between the observation and a random sample from the likelihood. The first argument must be the parameter vector.
  • prior: A Distribution defining the prior.
  • args...: Further positional arguments passed to f_dist
  • n_particles: Desired number of particles.
  • n_simulation: Maximal number of simulations from f_dist.
  • propsal = DifferentialEvolution(n_para = length(prior)): Method to generate propsals. Currently RandomWalk, DifferentialEvolution, and StretchMove are implemented.
  • algorithm = :single_eps: Algorithm for tolerance, either :multi_eps, or :single_eps. See below for details.
  • resample: After how many accepted population updates?
  • v = 1.0: Tuning parameter for annealing speed. Must be positive.
  • δ = 0.1: Tuning parameter for resampling intensity. Must be positive and should be small.
  • checkpoint_history = 1: every how many population updates distances and epsilons are stored
  • show_progressbar::Bool = !is_logging(stderr): defaults to true for interactive use.
  • show_checkpoint::Int = 100: every how many population updates algorithm state is displayed. By default disabled for for interactive use.
  • kwargs...: Further keyword arguments passed to f_dist`

Details

Depending on how many statistics f_dist returns, different algorithms are compatible:

1 statistic>1 statistics
:single_eps
:multi_eps

Note, there is no check if the chosen algorithm is compatible with f_dist!

Return

  • An object of type SABCresult
source
SimulatedAnnealingABC.update_population!Function
update_population!(population_state::SABCresult,
                   f_dist, prior, args...;
                   n_simulation,
                   v=1.0, δ=0.1,
                   proposal::Proposal = DifferentialEvolution(n_para = length(prior)),
                   resample = 2*length(population_state.population),
                   checkpoint_history = 1,
                   show_progressbar::Bool = !is_logging(stderr),
                   show_checkpoint = is_logging(stderr) ? 100 : Inf,
                   kwargs...)

Updates particles with n_simulation and applies importance sampling if needed. Modifies population_state.

Arguments

See docstring for sabc.

source
SimulatedAnnealingABC.SABCresultType

Holds results from a SABC run with fields:

  • population: vector of parameter samples from the approximate posterior
  • u: transformed distances
  • ρ: distances
  • state: state of algorithm

The history of ϵ can be accessed with the field state.ϵ_history. The history of ρ can be accessed with the field state.ρ_history. The history of u can be accessed with the field state.u_history.

source

Proposals

SimulatedAnnealingABC.DifferentialEvolutionType
DifferentialEvolution(; n_para, σ_gamma = 1e-5)
DifferentialEvolution(; γ0, σ_gamma = 1e-5)

Differential Evolution proposal, default values corresponding to EMCEE. If the number of parameters n_para is provided, γ0 is set to 2.38 / sqrt(2 * n_parameters).

References

Ter Braak, C.J., 2006. A Markov Chain Monte Carlo version of the genetic algorithm Differential Evolution: easy Bayesian computing for real parameter spaces. Statistics and Computing 16, 239–249.

Nelson, B., Ford, E.B., Payne, M.J., 2013. Run Dmc: An Efficient, Parallel Code For Analyzing Radial Velocity Observations Using N-Body Integrations And Differential Evolution Markov Chain Monte Carlo. ApJS 210, 11. https://doi.org/10.1088/0067-0049/210/1/11

source
SimulatedAnnealingABC.StretchMoveType
StretchMove(;a=2)

The standard proposal used in EMCEE.

Reference

Goodman, J., Weare, J., 2010. Ensemble samplers with affine invariance. Communications in Applied Mathematics and Computational Science 5, 65–80.

source
SimulatedAnnealingABC.RandomWalkType
RandomWalk(; β=0.8, n_para)

Gaussian random walk proposal.

The covariance is adaptivily learned. The mixing is controlled by the tuning parameter β which must be between zero and one.

source