Tempering MCMC

Tempering MCMC algorithms aim at sampling from a target distribution \(p_1(x)\) of the form \(p_1(x)=\frac{q_1(x)p_0(x)}{Z_1}\) where the factor \(q_1(x)\) and reference distribution \(p_0(x)\) can be evaluated. Additionally, these algorithms return an estimate of the normalization constant \(Z_1=\int{q_{1}(x) p_{0}(x)dx}\).

The algorithms sample from a sequence of intermediate densities \(p_{\beta}(x) \propto q_{\beta}(x) p_{0}(x)\) for values of the parameter \(\beta\) between 0 and 1 (\(\beta=\frac{1}{T}\) where \(T\) is sometimes called the temperature, \(q_{\beta}(x)\) is referred to as the intermediate factor associated with tempering parameter \(\beta\)). Setting \(\beta = 1\) equates sampling from the target, while \(\beta \rightarrow 0\) samples from the reference distribution.

Parallel tempering samples from all distributions simultaneously, and the tempering parameters \(0 < \beta_1 < \beta_2 < \cdots < \beta_{N} \leq 1\) must be chosen in advance by the user. Sequential tempering on the other hand samples from the various distributions sequentially, starting from the reference distribution, and the tempering parameters are selected adaptively by the algorithm.

The TemperingMCMC base class defines inputs that are common to parallel and sequential tempering:

class TemperingMCMC(pdf_intermediate=None, log_pdf_intermediate=None, args_pdf_intermediate=(), distribution_reference=None, save_log_pdf=True, random_state=None)[source]

Parent class to parallel and sequential tempering MCMC algorithms.

Parameters:
  • pdf_intermediate – callable that computes the intermediate factor. It should take at least two inputs x (ndarray, point(s) at which to evaluate the function), and temper_param (float, tempering parameter). Eit her pdf_intermediate or log_pdf_intermediate must be provided (log_pdf_intermediate is preferred). Within the code, the log_pdf_intermediate is evaluated as: log_pdf_intermediate(x, temper_param, *args_pdf_intermediate) where args_pdf_intermediate are additional positional arguments that are provided to the class via its args_pdf_intermediate input

  • log_pdf_intermediate – see pdf_intermediate

  • args_pdf_intermediate – see pdf_intermediate

  • distribution_reference – reference pdf \(p_0\) as a Distribution object

  • save_log_pdf – see same input in MCMC

Parallel Tempering

This algorithm (see e.g. [33] for theory about parallel tempering) runs the chains sampling from the various tempered distributions simultaneously. Periodically during the run, the different temperatures swap members of their ensemble in a way that preserves detailed balance. The chains closer to the reference chain (hot chains) can sample from regions that have low probability under the target and thus allow a better exploration of the parameter space, while the cold chains can better explore regions of high likelihood.

In parallel tempering, the normalizing constant \(Z_1\) is evaluated via thermodynamic integration ([34]). Defining the potential function \(U_{\beta}(x)=\frac{\partial \log{q_{\beta}(x)}}{\partial \beta}\), then

\(\ln{Z_1} = \log{Z_0} + \int_{0}^{1} E_{x \sim p_{\beta}} \left[ U_{\beta}(x) \right] d\beta\)

where the expectations are approximated via MC sampling using saved samples from the intermediate distributions. A trapezoidal rule is used for integration.

The ParallelTemperingMCMC class is imported using the following command:

>>> from UQpy.sampling.mcmc.tempering_mcmc.ParallelTemperingMCMC import ParallelTemperingMCMC
class ParallelTemperingMCMC(n_iterations_between_sweeps, pdf_intermediate=None, log_pdf_intermediate=None, args_pdf_intermediate=(), distribution_reference=None, save_log_pdf=False, nsamples=None, nsamples_per_chain=None, random_state=None, tempering_parameters=None, n_tempering_parameters=None, samplers=None)[source]

Class for Parallel-Tempering MCMC.

Parameters:
  • save_log_pdf (bool) – boolean, see MCMC documentation. Importantly, this needs to be set to True if one wants to evaluate the normalization constant via thermodynamic integration.

  • n_iterations_between_sweeps (int) – number of iterations (sampling steps) between sweeps between chains.

  • tempering_parameters (Optional[list]) – tempering parameters, as a list of N floats increasing from 0. to 1. Either tempering_parameters or n_tempering_parameters should be provided

  • n_tempering_parameters (Optional[int]) – number of tempering levels N, the tempering parameters are selected to follow a geometric suite by default

  • samplers (Union[MCMC, list[MCMC], None]) – MCMC object or list of such objects: MCMC samplers used to sample the parallel chains. If only one object is provided, the same MCMC sampler is used for all chains. Default to running a simple MH algorithm, where the proposal covariance for a given chain is inversely proportional to the tempering parameter.

run(nsamples=None, nsamples_per_chain=None)[source]

Run the MCMC algorithm.

This function samples from the MCMC chains and appends samples to existing ones (if any). This method leverages the run_iterations method specific to each of the samplers.

Parameters:
  • nsamples (Optional[int]) – Number of samples to generate from the target (the same number of samples will be generated for all intermediate distributions).

  • nsamples_per_chain (Optional[int]) – Number of samples per chain to generate from the target. Either nsamples or nsamples_per_chain must be provided (not both)

evaluate_normalization_constant(compute_potential, log_Z0=None, nsamples_from_p0=None)[source]

Evaluate normalization constant \(Z_1\).

The function returns an approximation of \(Z_1\), and saves intermediate results (value of \(\ln{Z_0}\), list of tempering parameters used in integration, and values of the associated expected potentials.

Parameters:
  • compute_potential – Function that takes three inputs: x (sample points where to evaluate the potential), log_factor_tempered_values (values of the log intermediate factors evaluated at points x), temper_param (tempering parameter) and evaluates the potential:

  • log_Z0 (Optional[float]) – Value of \(\ln{Z_{0}}\) (float), if unknwon, see nsamples_from_p0.

  • nsamples_from_p0 (Optional[int]) – number of samples from the reference distribution to sample to evaluate \(\ln{Z_{0}}\). Used only if input log_Z0 is not provided.

Sequential Tempering

This algorithm (first introduced in [35]) samples from a series of intermediate targets that are each tempered versions of the final/true target. In going from one intermediate distribution to the next, the existing samples are resampled according to some weights (similar to importance sampling). To ensure that there aren’t a large number of duplicates, the resampling step is followed by a short (or even single-step) Metropolis Hastings run that disperses the samples while remaining within the correct intermediate distribution. The final intermediate target is the required target distribution, and the samples following this distribution are the required samples.

The normalization constant \(Z_1\) is estimated as the product of the normalized sums of the resampling weights for each intermediate distribution, i.e. if \(w_{\beta_j}(x_{j_i})\) is the resampling weight corresponding to tempering parameter \(\beta_j\), calculated for the i-th sample for the intermediate distribution associated with \(\beta_j\), then \(Z_1 = \prod_{j=1}^{N} \ [ \sum_{i=i}^{\text{nsamples}} \ ]\). The Coefficient of Variance (COV) for this estimator is also given in [35].

The SequentialTemperingMCMC class is imported using the following command:

>>> from UQpy.sampling.mcmc.tempering_mcmc.SequentialTemperingMCMC import SequentialTemperingMCMC
class SequentialTemperingMCMC(pdf_intermediate=None, log_pdf_intermediate=None, args_pdf_intermediate=(), seed=None, distribution_reference=None, sampler=None, nsamples=None, recalculate_weights=False, save_intermediate_samples=False, percentage_resampling=100, random_state=None, resampling_burn_length=0, resampling_proposal=None, resampling_proposal_is_symmetric=True)[source]

Class for Sequential-Tempering MCMC

Parameters:
  • sampler (Optional[MCMC]) – MCMC object: MCMC samplers used to draw the remaining samples for the intermediate distribution after the resampling step. Default to running a simple MH algorithm, where the proposal covariance is calculated as per the procedure given in Ching and Chen (2007) and Betz et. al. (2016).

  • recalculate_weights (bool) – boolean: To be set to true if the resampling weights are to be recalculated after each point is generated during the resampling step. This is done so that the resampling weights are in accordance with the new sample generated after Metropolis Hastings is used for dispersion to ensure uniqueness of the samples.

  • save_intermediate_samples – boolean: To be set to true to save the samples that are generated according to the intermediate distributions.

  • percentage_resampling (int) – float: Indicates what percentage of samples for a given intermediate distribution are to be generated through resampling from the set of samples generated for the previous intermediate distribution.

  • resampling_burn_length (int) – int: Burn-in length for the Metropolis Hastings dispersion step to ensure uniqueness.

  • resampling_proposal (Optional[Distribution]) – Distribution object. The proposal distribution for the Metropolis Hastings dispersion step.

  • resampling_proposal_is_symmetric (bool) – boolean: Indicates whether the provided resampling proposal is symmetric.

run(nsamples=None)[source]

Run the MCMC algorithm.

This function samples from each intermediate distribution until samples from the target are generated. Samples cannot be appended to existing samples in this method. It leverages the run_iterations method specific to the sampler.

Parameters:

nsamples (Optional[int]) – Number of samples to generate from the target (the same number of samples will be generated for all intermediate distributions).

evaluate_normalization_constant()[source]

Computes the normalization constant \(Z_{1}=\int{q_{1}(x) p_{0}(x)dx}\) where \(p_0\) is the reference pdf and \(q_1\) is the target factor.

Examples