Obsidian
Probabilistic Geophysical Joint Inversion
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
stateline::mcmc Namespace Reference

Classes

class  ChainArray
 Manager for all the states and handle reading / writing from / to databse. More...
 
class  EpsrConvergenceCriteria
 Convergence test using estimated potential scale reduction (EPSR). It is a convergence metric that takes into account the variance of the means between chains and the variance of the samples within each chain. More...
 
class  Sampler
 Represents a Markov Chain Monte Carlo sampler that returns samples from a distribution. More...
 
struct  State
 Used to represent a state in a Markov Chain. More...
 

Typedefs

using PropAcceptFn = std::function< bool(const State &, const State &, double)>
 Type representing proposal acceptance functions.
 
using SwapAcceptFn = std::function< bool(const State &, const State &, double, double)>
 Type representing swap acceptance functions.
 

Enumerations

enum  SwapType { SwapType::NoAttempt, SwapType::Accept, SwapType::Reject }
 Used for recording swapping of chains. More...
 

Functions

Eigen::VectorXd bouncyBounds (const Eigen::VectorXd &val, const Eigen::VectorXd &min, const Eigen::VectorXd &max)
 A function to bounce the MCMC proposal off the hard boundaries. This allows the proposal to always move around without getting stuck at 'walls'. More...
 
Eigen::VectorXd adaptiveGaussianProposal (const Eigen::VectorXd &state, double sigma, const Eigen::VectorXd &min, const Eigen::VectorXd &max)
 An adaptive Gaussian proposal function. It randomly varies each value in the state according to a Gaussian distribution whose variance changes depending on the acceptance ratio of a chain. It also bounces of the walls of the hard boundaries given so as not to get stuck in corners. More...
 
bool acceptProposal (const State &newState, const State &oldState, double beta)
 Returns true if we want to accept the MCMC step. More...
 
bool acceptSwap (const State &stateLow, const State &stateHigh, double betaLow, double betaHigh)
 Returns true if we want to accept the MCMC swap. More...
 

Detailed Description

Namespace for MCMC sampling features.

Enumeration Type Documentation

Used for recording swapping of chains.

Enumerator
NoAttempt 

No swap was attempted.

Accept 

A swap was accepted.

Reject 

A swap was rejected.

Function Documentation

bool stateline::mcmc::acceptProposal ( const State &  newState,
const State &  oldState,
double  beta 
)

Returns true if we want to accept the MCMC step.

Parameters
newStateThe proposed state.
oldStateThe current state of the chain.
betaThe inverse temperature of the chain.
Returns
True if the proposal was accepted.
bool stateline::mcmc::acceptSwap ( const State &  stateLow,
const State &  stateHigh,
double  betaLow,
double  betaHigh 
)

Returns true if we want to accept the MCMC swap.

Parameters
stateLowThe state of the lower temperature chain.
stateHighThe state of the higher temperature chain.
betaLowThe inverse temperature of the lower temperature chain.
betaHighThe inverse temperature of the high temperature chain.
Returns
True if the swap was accepted.
Eigen::VectorXd stateline::mcmc::adaptiveGaussianProposal ( const Eigen::VectorXd &  state,
double  sigma,
const Eigen::VectorXd &  min,
const Eigen::VectorXd &  max 
)

An adaptive Gaussian proposal function. It randomly varies each value in the state according to a Gaussian distribution whose variance changes depending on the acceptance ratio of a chain. It also bounces of the walls of the hard boundaries given so as not to get stuck in corners.

Parameters
stateThe current state of the chain
sigmaThe standard deviation of the distribution (step size of the proposal)
minThe minimum bound of theta
maxThe maximum bound of theta
Returns
The new proposed theta
Eigen::VectorXd stateline::mcmc::bouncyBounds ( const Eigen::VectorXd &  val,
const Eigen::VectorXd &  min,
const Eigen::VectorXd &  max 
)

A function to bounce the MCMC proposal off the hard boundaries. This allows the proposal to always move around without getting stuck at 'walls'.

Parameters
valThe proposed value
minThe minimum bound of theta
maxThe maximum bound of theta
Returns
The new bounced theta definitely in the bounds