API/Reference
NestedSamplers.Bounds
NestedSamplers.Models
NestedSamplers.Proposals
NestedSamplers.Bounds.Ellipsoid
NestedSamplers.Bounds.MultiEllipsoid
NestedSamplers.Bounds.NoBounds
NestedSamplers.Nested
NestedSamplers.NestedModel
NestedSamplers.Proposals.RSlice
NestedSamplers.Proposals.RStagger
NestedSamplers.Proposals.RWalk
NestedSamplers.Proposals.Slice
NestedSamplers.Proposals.Uniform
NestedSamplers.Models.CorrelatedGaussian
NestedSamplers.Models.GaussianShells
Samplers
NestedSamplers.NestedModel
— TypeNestedModel(loglike, prior_transform)
NestedModel(loglike, priors::AbstractVector{<:Distribution})
loglike
must be callable with a signature loglike(::AbstractVector)
where the length of the vector must match the number of parameters in your model.
prior_transform
must be a callable with a signature prior_transform(::AbstractVector)
that returns the transformation from the unit-cube to parameter space. This is effectively the quantile or ppf of a statistical distribution. For convenience, if a vector of Distribution
is provided (as a set of priors), a transformation function will automatically be constructed using Distributions.quantile
.
Note: loglike
is the only function used for likelihood calculations. This means if you want your priors to be used for the likelihood calculations they must be manually included in the loglike
function.
NestedSamplers.Nested
— TypeNested(ndims, nactive;
bounds=Bounds.MultiEllipsoid,
proposal=:auto,
enlarge=1.25,
update_interval=default_update_interval(proposal, ndims),
min_ncall=2nactive,
min_eff=0.10)
Static nested sampler with nactive
active points and ndims
parameters.
ndims
is equivalent to the number of parameters to fit, which defines the dimensionality of the prior volume used in evidence sampling. nactive
is the number of live or active points in the prior volume. This is a static sampler, so the number of live points will be constant for all of the sampling.
Bounds and Proposals
bounds
declares the Type of Bounds.AbstractBoundingSpace
to use in the prior volume. The available bounds are described by Bounds
. proposal
declares the algorithm used for proposing new points. The available proposals are described in Proposals
. If proposal
is :auto
, will choose the proposal based on ndims
ndims < 10
-Proposals.Uniform
10 ≤ ndims ≤ 20
-Proposals.RWalk
ndims > 20
-Proposals.Slice
The original nested sampling algorithm is roughly equivalent to using Bounds.Ellipsoid
with Proposals.Uniform
. The MultiNest algorithm is roughly equivalent to Bounds.MultiEllipsoid
with Proposals.Uniform
. The PolyChord algorithm is roughly equivalent to using Proposals.RSlice
.
Other Parameters
enlarge
- When fitting the bounds to live points, they will be enlarged (in terms of volume) by this linear factor.update_interval
- How often to refit the live points with the bounds as a fraction ofnactive
. By default this will be determined usingdefault_update_interval
for the given proposalProposals.Uniform
-1.5
Proposals.RWalk
andProposals.RStagger
-0.15 * walks
Proposals.Slice
-0.9 * ndims * slices
Proposals.RSlice
-2 * slices
min_ncall
- The minimum number of iterations before trying to fit the first boundmin_eff
- The maximum efficiency before trying to fit the first bound
Convergence
There are a few convergence criteria available, by default the dlogz
criterion will be used.
dlogz=0.5
sample until the fraction of the remaining evidence is below the given value (more info).maxiter=Inf
stop after the given number of iterationsmaxcall=Inf
stop after the given number of log-likelihood function callsmaxlogl=Inf
stop after reaching the target log-likelihood
Bounds
NestedSamplers.Bounds
— ModuleNestedSamplers.Bounds
This module contains the different algorithms for bounding the prior volume.
The available implementations are
Bounds.NoBounds
- no bounds on the prior volume (equivalent to a unit cube)Bounds.Ellipsoid
- bound using a single ellipsoidBounds.MultiEllipsoid
- bound using multiple ellipsoids in an optimal cluster
NestedSamplers.Bounds.NoBounds
— TypeBounds.NoBounds([T=Float64], N)
Unbounded prior volume; equivalent to the unit cube in N
dimensions.
NestedSamplers.Bounds.Ellipsoid
— TypeBounds.Ellipsoid([T=Float64], N)
Bounds.Ellipsoid(center::AbstractVector, A::AbstractMatrix)
An N
-dimensional ellipsoid defined by
\[(x - center)^T A (x - center) = 1\]
where size(center) == (N,)
and size(A) == (N,N)
.
NestedSamplers.Bounds.MultiEllipsoid
— TypeBounds.MultiEllipsoid([T=Float64], ndims)
Bounds.MultiEllipsoid(::AbstractVector{Ellipsoid})
Use multiple Ellipsoid
s in an optimal clustering to bound prior space. For more details about the bounding algorithm, see the extended help (??Bounds.MultiEllipsoid
)
Proposals
NestedSamplers.Proposals
— ModuleNestedSamplers.Proposals
This module contains the different algorithms for proposing new points within a bounding volume in unit space.
The available implementations are
Proposals.Uniform
- samples uniformly within the bounding volumeProposals.RWalk
- random walks to a new point given an existing oneProposals.RStagger
- random staggering away to a new point given an existing oneProposals.Slice
- slicing away to a new point given an existing oneProposals.RSlice
- random slicing away to a new point given an existing one
NestedSamplers.Proposals.Uniform
— TypeProposals.Uniform()
Propose a new live point by uniformly sampling within the bounding volume.
NestedSamplers.Proposals.RWalk
— TypeProposals.RWalk(;ratio=0.5, walks=25, scale=1)
Propose a new live point by random walking away from an existing live point.
Parameters
ratio
is the target acceptance ratiowalks
is the minimum number of steps to takescale
is the proposal distribution scale, which will update between proposals.
NestedSamplers.Proposals.RStagger
— TypeProposals.RStagger(;ratio=0.5, walks=25, scale=1)
Propose a new live point by random staggering away from an existing live point. This differs from the random walk proposal in that the step size here is exponentially adjusted to reach a target acceptance rate during each proposal, in addition to between proposals.
Parameters
ratio
is the target acceptance ratiowalks
is the minimum number of steps to takescale
is the proposal distribution scale, which will update between proposals.
NestedSamplers.Proposals.Slice
— TypeProposals.Slice(;slices=5, scale=1)
Propose a new live point by a series of random slices away from an existing live point. This is a standard Gibbs-like implementation where a single multivariate slice is a combination of slices
univariate slices through each axis.
Parameters
slices
is the minimum number of slicesscale
is the proposal distribution scale, which will update between proposals.
NestedSamplers.Proposals.RSlice
— TypeProposals.RSlice(;slices=5, scale=1)
Propose a new live point by a series of random slices away from an existing live point. This is a standard random implementation where each slice is along a random direction based on the provided axes.
Parameters
slices
is the minimum number of slicesscale
is the proposal distribution scale, which will update between proposals.
Models
NestedSamplers.Models
— ModuleThis module contains various statistical models in the form of NestedModel
s. These models can be used for examples and for testing.
NestedSamplers.Models.GaussianShells
— FunctionModels.GaussianShells()
2-D Gaussian shells centered at [-3.5, 0]
and [3.5, 0]
with a radius of 2 and a shell width of 0.1
Examples
julia> model, lnZ = Models.GaussianShells();
julia> lnZ
-1.75
NestedSamplers.Models.CorrelatedGaussian
— FunctionModels.CorrelatedGaussian(ndims)
Creates a highly-correlated Gaussian with the given dimensionality.
\[\mathbf\theta \sim \mathcal{N}\left(2\mathbf{1}, \mathbf{I}\right)\]
\[\Sigma_{ij} = \begin{cases} 1 &\quad i=j \\ 0.95 &\quad i\neq j \end{cases}\]
\[\mathcal{L}(\mathbf\theta) = \mathcal{N}\left(\mathbf\theta | \mathbf{0}, \mathbf\Sigma \right)\]
the analytical evidence of the model is
\[Z = \mathcal{N}\left(2\mathbf{1} | \mathbf{0}, \mathbf\Sigma + \mathbf{I} \right)\]
Examples
julia> model, lnZ = Models.CorrelatedGaussian(10);
julia> lnZ
-12.482738597926607