Nested rejection sampler

class pints.NestedRejectionSampler(log_prior)[source]

Creates a nested sampler that estimates the marginal likelihood and generates samples from the posterior.

This is the simplest form of nested sampler and involves using rejection sampling from the prior as described in the algorithm on page 839 in [1] to estimate the marginal likelihood and generate weights, preliminary samples (with their respective likelihoods), required to generate posterior samples.

The posterior samples are generated as described in [1] on page 849 by randomly sampling the preliminary point, accounting for their weights and likelihoods.

Initialise:

Z = 0
X_0 = 1

Draw samples from prior:

for i in 1:n_active_points:
    theta_i ~ p(theta), i.e. sample from the prior
    L_i = p(theta_i|X)
endfor

In each iteration of the algorithm (t):

L_min = min(L)
indexmin = min_index(L)
X_t = exp(-t / n_active_points)
w_t = X_t - X_t-1
Z = Z + L_min * w_t
theta* ~ p(theta)
while p(theta*|X) < L_min:
    theta* ~ p(theta)
endwhile
theta_indexmin = theta*
L_indexmin = p(theta*|X)

At the end of iterations, there is a final Z increment:

Z = Z + (1 / n_active_points) * (L_1 + L_2 + ..., + L_n_active_points)

The posterior samples are generated as described in [1] on page 849 by weighting each dropped sample in proportion to the volume of the posterior region it was sampled from. That is, the probability for drawing a given sample j is given by:

p_j = L_j * w_j / Z

where j = 1, …, n_iterations.

Extends NestedSampler.

References

[1](1, 2) “Nested Sampling for General Bayesian Computation”, John Skilling, Bayesian Analysis 1:4 (2006). https://doi.org/10.1214/06-BA127
active_points()

Returns the active points from nested sampling run.

ask(n_points)[source]

Proposes new point(s) by sampling from the prior.

in_initial_phase()

For methods that need an initial phase (see needs_initial_phase()), this method returns True if the method is currently configured to be in its initial phase. For other methods a NotImplementedError is returned.

min_index()

Returns index of sample with lowest log-likelihood.

n_active_points()

Returns the number of active points that will be used in next run.

n_hyper_parameters()[source]

See TunableMethod.n_hyper_parameters().

name()[source]

See pints.NestedSampler.name().

needs_initial_phase()

Returns True if this method needs an initial phase, for example ellipsoidal nested sampling has a period of running rejection sampling before it starts to fit ellipsoids to points.

needs_sensitivities()

Determines whether sampler uses sensitivities of the solution.

running_log_likelihood()

Returns current value of the threshold log-likelihood value.

set_hyper_parameters(x)[source]

Hyper-parameter vector is: [active_points_rate]

Parameters:x – An array of length n_hyper_parameters used to set the hyper-parameters
set_initial_phase(in_initial_phase)

For methods that need an initial phase (see needs_initial_phase()), this method toggles the initial phase algorithm. For other methods a NotImplementedError is returned.

set_n_active_points(active_points)

Sets the number of active points for the next run.

tell(fx)

If a single evaluation is provided as arguments, a single point is accepted and returned if its likelihood exceeds the current threshold; otherwise None is returned.

If multiple evaluations are provided as arguments (for example, if running the algorithm in parallel), None is returned if no points have likelihood exceeding threshold; if a single point passes the threshold, it is returned; if multiple points pass, one is selected uniformly at random and returned and the others are stored for later use.

In all cases, two objects are returned: the proposed point (which may be None) and an array of other points that also pass the threshold (which is empty for single evaluation mode but may be non-empty for multiple evaluation mode).