mcmc

mcmc

Base module for MCMC methods.

class quinn.mcmc.mcmc.MCMCBase[source]

Bases: object

Base class for MCMC.

logPost

Log-posterior evaluator function. It has a singature of logPost(model_parameters, **postInfo) returning a float.

Type:

callable

logPostGrad

Log-posterior gradient evaluator function. It has a singature of logPostGrad(model_parameters, **postInfo) returning an np.ndarray of size model_parameters. Defaults to None, for non-gradient based methods.

Type:

callable

postInfo

Dictionary that holds auxiliary parameters for posterior evaluations.

Type:

dict

__init__()[source]

Dummy instantiation.

setLogPost(logPost, logPostGrad, **postInfo)[source]

Setting logposterior and its gradient functions.

Parameters:
  • logPost (callable) – Log-posterior evaluator function. It has a singature of logPost(model_parameters, **postInfo) returning a float.

  • logPostGrad (callable) – Log-posterior gradient evaluator function. It has a singature of logPostGrad(model_parameters, **postInfo) returning an np.ndarray of size model_parameters. Can be None, for non-gradient based methods.

  • postInfo (dict) – Dictionary that holds auxiliary parameters for posterior evaluations.

run(nmcmc, param_ini)[source]

Markov chain Monte Carlo running function.

Parameters:
  • nmcmc (int) – Number of steps.

  • param_ini (np.ndarray) – Initial state of the chain.

Returns:

Dictionary of results. Keys are ‘chain’ (chain samples array), ‘mapparams’ (MAP parameters array), ‘maxpost’ (maximal log-post value), ‘accrate’ (acceptance rate), ‘logpost’ (log-posterior array), ‘alphas’ (array of Metropolis-Hastings probability ratios).

Return type:

dict

sampler(current, imcmc)[source]

Sampler method.

Parameters:
  • current (np.ndarray) – Current chain state.

  • imcmc (int) – Current chain step number.

Raises:

NotImplementedError – Not implemented in the base class.

admcmc

Module for Adaptive MCMC (AMCMC) sampling. For details of the method, see Haario et al. [3]

class quinn.mcmc.admcmc.AMCMC(cov_ini=None, gamma=0.1, t0=100, tadapt=1000)[source]

Bases: MCMCBase

Adaptive MCMC class.

cov_ini

Initial covariance array of size (p,p).

Type:

np.ndarray

gamma

Proposal jump size factor \(\gamma\).

Type:

float

_propcov

A 2d array of size (p,p) for proposal covariance.

Type:

np.ndarray

t0

Step where adaptivity begins.

Type:

int

tadapt

Frequency for adapting/updating the covariance.

Type:

int

__init__(cov_ini=None, gamma=0.1, t0=100, tadapt=1000)[source]

Initialization.

Parameters:
  • cov_ini (np.ndarray, optional) – Initial covariance array of size (p,p). Defaults to None which sets the initial covariance as some fraction of the chain state.

  • gamma (float, optional) – Proposal jump size factor \(\gamma\). Defaults to None.

  • t0 (int, optional) – Step where adaptivity begins. Defaults to 100.

  • tadapt (int) – Frequency for adapting/updating the covariance. Defaults to 1000.

sampler(current, imcmc)[source]

Sampler method.

Parameters:
  • current (np.ndarray) – Current chain state.

  • imcmc (int) – Current chain step number.

Returns:

A tuple containing the current proposal sample, and two zeros irrelevant for AMCMC.

Return type:

tuple(current_proposal, current_K, proposed_K)

hmc

Module for Hamiltonian MCMC (HMC) sampling. For details of the method, see Chapter 5 of Brooks et al. [4] or https://arxiv.org/pdf/1206.1901.pdf

class quinn.mcmc.hmc.HMC(epsilon=0.05, L=3)[source]

Bases: MCMCBase

Hamiltonian MCMC class.

epsilon

Step size of the method.

Type:

float

L

Number of steps in the Hamiltonian integrator.

Type:

int

__init__(epsilon=0.05, L=3)[source]

Initialization.

Parameters:
  • epsilon (float, optional) – Step size of the method. Defaults to 0.05.

  • L (int, optional) – Number of steps in the Hamiltonian integrator. Defaults to 3.

sampler(current, imcmc)[source]

Sampler method.

Parameters:
  • current (np.ndarray) – Current chain state.

  • imcmc (int) – Current chain step number.

Returns:

A tuple containing the current proposal sample, current and proposed kinetic energies.

Return type:

tuple(current_proposal, current_K, proposed_K)

mala

Module for Metropolis Adjusted Langevin (MALA) MCMC sampling. For details of the method, see Girolami and Calderhead [5] or https://arxiv.org/pdf/1206.1901.pdf

class quinn.mcmc.mala.MALA(epsilon=0.05)[source]

Bases: MCMCBase

MALA MCMC class.

epsilon

Step size of the method.

Type:

float

__init__(epsilon=0.05)[source]

Initialization.

Parameters:

epsilon (float, optional) – Step size of the method. Defaults to 0.05.

sampler(current, imcmc)[source]

Sampler method.

Parameters:
  • current (np.ndarray) – Current chain state.

  • imcmc (int) – Current chain step number.

Returns:

A tuple containing the current proposal sample, current and proposed kinetic energies.

Return type:

tuple(current_proposal, current_K, proposed_K)

Note: When the dust settles, MALA is actually exactly HMC with L=1.