MCMC
mcmc
¶
MLDA-specific MCMC routines.
This module contains all MCMC-specific functionalities, including those necessary in the Multilevel context. Roughly speaking, this comprises proposals for coarse level changes and accept-reject routines. The only prefetching-specific component is the accept rate estimator, which is used to predict possible future states of the Markov chain. Proposals and accept rate estimators are implemented in class hierarchies. New options can easily be implemented by inheriting from the respective base class interfaces.
Classes:
Name | Description |
---|---|
BaseProposal |
Base class for MCMC proposals |
RandomWalkProposal |
Metropolis random walk proposal |
PCNProposal |
Preconditioned Crank-Nicolson proposal |
BaseAcceptRateEstimator |
Base class for MLDA accept rate estimators |
StaticAcceptRateEstimator |
Static accept rate estimator |
MLMetropolisHastingsKernel |
Metropolis-Hastings acceptance kernel for multilevel decisions |
mtmlda.core.mcmc.BaseProposal
¶
Base class for MCMC proposals.
This class defines the interface for MCMC proposals. It is an abstract class and cannot be instantiated. It provides two abstract methods that need to be implemented by subclasses: propose and evaluate_log_probability. The propose method generates a new proposal based on the current state. The evaluate_log_probability method evaluates the log probability of a new move depending on the current state.
Methods:
Name | Description |
---|---|
propose |
Propose new MCMC move |
evaluate_log_probability |
Evaluate log probability of a new move depending on the current state |
Attributes: rng: Random number generator used for proposals
rng
property
writable
¶
Getter of the random number generator for proposals.
Returns:
Type | Description |
---|---|
np.random.Generator
|
np.random.Generator: Numpy RNG object |
__init__
¶
Base class constructor.
The constructor simply initializes a numpy random number generator using the provided seed.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
seed
|
int
|
Seed value for the random number generator |
required |
Raises:
Type | Description |
---|---|
TypeError
|
Checks for valid seed type |
propose
abstractmethod
¶
Abstract method for proposal generation.
This method needs to be re-implemented in child classes.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
current_state
|
np.ndarray
|
State of the Markov chain to propose from |
required |
evaluate_log_probability
abstractmethod
¶
Abstract method for log probability evaluation.
This method needs to be re-implemented in child classes. The states are not called "current" and "new" state, as their relation is not as clear in the multilevel setting. The names simply refer to the order of arguments in the evaluation formula.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
left_state
|
np.ndarray
|
First state for evaluation |
required |
right_state
|
np.ndarray
|
Second state for evaluation |
required |
mtmlda.core.mcmc.RandomWalkProposal
¶
Bases: BaseProposal
Random walk proposal.
Implementation of the Metropolis random walk proposal, inheriting from the BaseProposal class.
The proposal distribution is a Gaussian about the current state. We utilize a fixed covariance
matrix and step width. Note that the proposal distribution is symmetric w.r.t. to its arguments,
so that it vanishes in the Metropolis-Hastings acceptance ratio. Therefore, the
evaluate_log_probability
method only returns 0 as a dummy value.
__init__
¶
Proposal constructor.
The constructor takes a fixed step width and covariance matrix for the Gaussian proposal.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
step_width
|
float
|
Proposal step width , needs to be positive |
required |
covariance
|
np.ndarray
|
Proposal covariance matrix |
required |
seed
|
int
|
Seed for initialization of the proposal RNG |
required |
Raises:
Type | Description |
---|---|
ValueError
|
Checks if step width is a positive real number |
ValueError
|
Checks if covariance matrix can be Cholesky factorized |
propose
¶
Propose new move from current state.
Given the covariance :math:\\Sigma
and current state :math:m
, the new state is
proposed from the distribution :math:\\mathcal{N}(m, \\Sigma)
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
current_state
|
np.ndarray
|
Current state of the Markov chain to propose from |
required |
Returns:
Type | Description |
---|---|
np.ndarray
|
np.ndarray: Proposal state |
evaluate_log_probability
¶
Evaluate log probability of proposal, given current state.
This is only a dummy returning zero, as the conditional proposal probabilities vanish in the Metropolis-Hastings acceptance ratio.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
_proposal
|
np.ndarray
|
Proposal state |
required |
_current_state
|
np.ndarray
|
Current state |
required |
Returns:
Name | Type | Description |
---|---|---|
float |
float
|
Dummy log-probability, here always zero |
mtmlda.core.mcmc.PCNProposal
¶
Bases: BaseProposal
Preconditioned Crank-Nicolson random walk proposal.
Implementation of the pCN proposal, inheriting from the BaseProposal class. The pCN proposal is an asymmetric extension of the Metropolis random walk proposal. It has been developed particularly for high-dimensional parameter spaces.
__init__
¶
Proposal constructor.
The parameter beta
is analogous to the step with in the Metropolis random walk proposal.
It determines the "explicitness" of the underlying CN scheme.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
beta
|
float
|
Step width, in (0,1) |
required |
covariance
|
np.ndarray
|
Covariance matrix for Gaussian |
required |
seed
|
int
|
RNG seed |
required |
Raises:
Type | Description |
---|---|
ValueError
|
Checks that |
ValueError
|
Checks that covariance matrix can be Cholesky factorized |
ValueError
|
Checks that covariance matrix can be inverted |
propose
¶
Propose a new state from the current one.
Given the covariance :math:\\Sigma
and current state :math:m
, the new state is
proposed from the distribution :math:\\mathcal{N}(\\sqrt{1-\\beta^2}m, \\beta^2\\Sigma)
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
current_state
|
np.ndarray
|
Current state of the Markov chain |
required |
Returns:
Type | Description |
---|---|
np.ndarray
|
np.ndarray: New proposal |
evaluate_log_probability
¶
Evaluate the log probability of the proposal given the current state.
The log probability for a proposal :math:m'
from current state :math:m
is given by
.. math::
\Delta m = m' - \sqrt{1-\beta^2}m,
\log p = -\frac{1}{2} \Delta m^T \beta^2 \Sigma^{-1} \Delta m.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
proposal
|
np.ndarray
|
Proposal state |
required |
current_state
|
np.ndarray
|
Current state |
required |
Returns:
Name | Type | Description |
---|---|---|
float |
float
|
Log probability of the proposal given the current state |
mtmlda.core.mcmc.BaseAcceptRateEstimator
¶
Base class for MLDA accept rate estimators.
This class defines the interface for accept rate estimators. It is an abstract class and cannot be instantiated. The prefetching approach to our MLDA implementation requires a-priori estimates for accepting MCMC moves on the different levels of the model hierarchy. Such estimates can depend on any number of factors like the level itself, the current state, or the sample history. In the simplest case, the acceptance rates are constant.
Methods:
Name | Description |
---|---|
get_acceptance_rate |
Get the acceptance rate estimate for a given MCMC move |
get_acceptance_rate
abstractmethod
¶
Abstract getter for the acceptance rate.
mtmlda.core.mcmc.StaticAcceptRateEstimator
¶
Bases: BaseAcceptRateEstimator
Static accept rate estimator.
This class implements a simple accept rate estimator based on a fixed update scheme. The object is initialized with a list of estimates for each levels. The estimates are then incremented or decremented after each MCMC decision, depending on the acceptance of the move.
Methods:
Name | Description |
---|---|
get_acceptance_rate |
Get the acceptance rate estimate for a given MCMC move |
update |
Update the acceptance rate estimate after a MCMC move |
__init__
¶
Estimator Constructor.
The constructor takes a list of initial guesses for each level and a fixed parameter for the update rule.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
initial_guess
|
list[float]
|
List of initial guesses for the acceptance rates of each level |
required |
update_parameter
|
float
|
Update parameter |
required |
Raises:
Type | Description |
---|---|
ValueError
|
Checks that initial guesses are in (0,1) |
ValueError
|
CHecks that update parameter is in (0,1) |
get_acceptance_rate
¶
Return acceptance rate.
For this simple estimator, the acceptance rate only depends on the current level.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
level
|
int
|
Level in the MLDA algorithm |
required |
Returns:
Name | Type | Description |
---|---|---|
float |
float
|
Acceptance rate estimate |
update
¶
Update an acceptance rate estimate.
Given a node, the update modifies the estimate corresponding to the node's level.
Let :math:\\alpha
be the current estimate and :math:\\delta
the update parameter.
If the corresponding MCMC move has been accepted, the new estimate is
:math:(1-\\delta)\\alpha + \\delta
. Otherwise, the new estimate is
:math:(1-\\delta)\\alpha
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
accepted
|
bool
|
Bool telling if the MCMC move has been accepted |
required |
node
|
mltree.MTNode
|
Node from which the move has been performed |
required |
mtmlda.core.mcmc.MLMetropolisHastingsKernel
¶
Metropolis-Hastings Acceptance Kernel.
The kernel implements the Metropolis-Hastings acceptance rule for the multilevel setting. Accordingly, to different acceptance rules are implemented. One for within-level moves, utilized for standard MCMC on the coarsest level of the model hierarchy. The other is for between-level moves on the higher levels.
Methods:
Name | Description |
---|---|
compute_single_level_decision |
Compute the decision for a within-level MCMC move |
compute_two_level_decision |
Compute the decision for a between-level MCMC move |
__init__
¶
Kernel constructor.
Takes an object derived from BaseProposal.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
ground_proposal
|
BaseProposal
|
Proposal object |
required |
compute_single_level_decision
¶
compute_two_level_decision
staticmethod
¶
Compute between-level MCMC decision.
This decision rule is specific to the multilevel setting. Note that it does not utilize the proposal, but only the node containing the proposal and its same-level parent.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
node
|
mltree.MTNode
|
Node containing the proposal |
required |
same_level_parent
|
mltree.MTNode
|
Same level parent of that node |
required |
Returns:
Name | Type | Description |
---|---|---|
bool |
bool
|
If proposal has been accepted |