MUQ  0.4.3
muq::SamplingAlgorithms::GMHKernel Class Reference

A kernel for the generalized Metropolis-Hastings kernel. More...

#include <GMHKernel.h>

Inheritance diagram for muq::SamplingAlgorithms::GMHKernel:

Detailed Description

A kernel for the generalized Metropolis-Hastings kernel.

Reference: "A general construction for parallelizing Metropolis-Hastings algorithms" (Calderhead, 2014)

Definition at line 19 of file GMHKernel.h.

Public Member Functions

 GMHKernel (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem)
 
 GMHKernel (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem, std::shared_ptr< MCMCProposal > proposalIn)
 
virtual ~GMHKernel ()
 
virtual void PreStep (unsigned int const t, std::shared_ptr< SamplingState > state) override
 
virtual std::vector< std::shared_ptr< SamplingState > > Step (unsigned int const t, std::shared_ptr< SamplingState > state) override
 
std::vector< std::shared_ptr< SamplingState > > SampleStationary () const
 Sample the stationary distribution. More...
 
Eigen::VectorXd StationaryAcceptance () const
 Get the cumulative stationary acceptance probability. More...
 
std::shared_ptr< parcer::Communicator > GetCommunicator () const
 
- Public Member Functions inherited from muq::SamplingAlgorithms::MHKernel
 MHKernel (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem)
 
 MHKernel (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem, std::shared_ptr< MCMCProposal > proposalIn)
 
virtual ~MHKernel ()=default
 
virtual std::shared_ptr< MCMCProposalProposal ()
 
virtual void PostStep (unsigned int const t, std::vector< std::shared_ptr< SamplingState >> const &state) override
 Allow the kernel to adapt given a new state. More...
 
virtual void PrintStatus (std::string prefix) const override
 
virtual double AcceptanceRate () const
 
virtual void SetBlockInd (int newBlockInd) override
 
- Public Member Functions inherited from muq::SamplingAlgorithms::TransitionKernel
 TransitionKernel (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem)
 
virtual ~TransitionKernel ()=default
 
virtual void SetCommunicator (std::shared_ptr< parcer::Communicator > newcomm)
 
virtual void PrintStatus () const
 
virtual int GetBlockInd () const
 
virtual std::shared_ptr< AbstractSamplingProblem > const & Problem () const
 

Additional Inherited Members

- Public Types inherited from muq::SamplingAlgorithms::TransitionKernel
typedef boost::function< std::shared_ptr< TransitionKernel >boost::property_tree::ptree, std::shared_ptr< AbstractSamplingProblem >)> TransitionKernelConstructor
 
typedef std::map< std::string, TransitionKernelConstructorTransitionKernelMap
 
- Static Public Member Functions inherited from muq::SamplingAlgorithms::TransitionKernel
static std::shared_ptr< TransitionKernelConstruct (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem)
 Static constructor for the transition kernel. More...
 
static std::shared_ptr< TransitionKernelMapGetTransitionKernelMap ()
 

Constructor & Destructor Documentation

◆ GMHKernel() [1/2]

muq::SamplingAlgorithms::GMHKernel::GMHKernel ( boost::property_tree::ptree const &  pt,
std::shared_ptr< AbstractSamplingProblem problem 
)
Parameters
[in]ptOptions for this kenel and the standard Metropolis-Hastings kernel
[in]problemThe problem we want to sample

◆ GMHKernel() [2/2]

muq::SamplingAlgorithms::GMHKernel::GMHKernel ( boost::property_tree::ptree const &  pt,
std::shared_ptr< AbstractSamplingProblem problem,
std::shared_ptr< MCMCProposal proposalIn 
)
Parameters
[in]ptOptions for this kenel and the standard Metropolis-Hastings kernel
[in]problemThe problem we want to sample
[in]proposalInThe proposal for the MCMC chain

◆ ~GMHKernel()

GMHKernel::~GMHKernel ( )
virtual

Definition at line 56 of file GMHKernel.cpp.

Member Function Documentation

◆ AcceptanceDensity()

void GMHKernel::AcceptanceDensity ( Eigen::VectorXd &  R)
private

Compute the stationary transition density.

Parameters
[in]RThe log-target

Definition at line 140 of file GMHKernel.cpp.

References ComputeStationaryAcceptance(), nlohmann::detail::dtoa_impl::k, Np1, muq::SamplingAlgorithms::MHKernel::proposal, and proposedStates.

Referenced by SerialProposal().

◆ AcceptanceMatrix()

Eigen::MatrixXd GMHKernel::AcceptanceMatrix ( Eigen::VectorXd const &  R) const
private

Definition at line 154 of file GMHKernel.cpp.

References Np1.

Referenced by ComputeStationaryAcceptance().

◆ ComputeStationaryAcceptance()

void GMHKernel::ComputeStationaryAcceptance ( Eigen::VectorXd const &  R)
private

Compute the cumulative acceptance density.

Parameters
[in]RThe log-target plus the log-proposals

Definition at line 168 of file GMHKernel.cpp.

References AcceptanceMatrix(), Np1, and stationaryAcceptance.

Referenced by AcceptanceDensity().

◆ GetCommunicator()

std::shared_ptr< parcer::Communicator > GMHKernel::GetCommunicator ( ) const

Definition at line 228 of file GMHKernel.cpp.

References muq::SamplingAlgorithms::TransitionKernel::comm.

◆ ParallelProposal()

void GMHKernel::ParallelProposal ( unsigned int const  t,
std::shared_ptr< SamplingState state 
)
private

Propose \(N\) points in parallel and evaluate the log target.

Parameters
[in]tThe current step in the MCMC chain
[in]stateThe current point

Definition at line 84 of file GMHKernel.cpp.

References muq::SamplingAlgorithms::TransitionKernel::comm, and SerialProposal().

Referenced by PreStep().

◆ PreStep()

void GMHKernel::PreStep ( unsigned int const  t,
std::shared_ptr< SamplingState state 
)
overridevirtual

Propose GMHKernel::N points and compute the cumulative distribution of the stationary distribution for the acceptance probability (GMHKernel::proposedStates and GMHKernel::stationaryAcceptance, respectively)

Parameters
[in]tThe current step in the MCMC chain
[in]stateThe current MCMC state

Reimplemented from muq::SamplingAlgorithms::TransitionKernel.

Definition at line 183 of file GMHKernel.cpp.

References muq::SamplingAlgorithms::TransitionKernel::comm, ParallelProposal(), and SerialProposal().

◆ SampleStationary()

std::vector< std::shared_ptr< SamplingState > > GMHKernel::SampleStationary ( ) const

Sample the stationary distribution.

Returns
The new points

Definition at line 197 of file GMHKernel.cpp.

References M, proposedStates, and stationaryAcceptance.

Referenced by Step().

◆ SerialProposal()

void GMHKernel::SerialProposal ( unsigned int const  t,
std::shared_ptr< SamplingState state 
)
private

Propose \(N\) points in serial and evaluate the log target.

Parameters
[in]tThe current step in the MCMC chain
[in]stateThe current point

Definition at line 58 of file GMHKernel.cpp.

References AcceptanceDensity(), Np1, muq::SamplingAlgorithms::TransitionKernel::problem, muq::SamplingAlgorithms::MHKernel::proposal, and proposedStates.

Referenced by ParallelProposal(), and PreStep().

◆ StationaryAcceptance()

Eigen::VectorXd GMHKernel::StationaryAcceptance ( ) const

Get the cumulative stationary acceptance probability.

Must be called after GMHKernel::PreStep.

Returns
The stationary acceptance probability

Definition at line 220 of file GMHKernel.cpp.

References Np1, and stationaryAcceptance.

◆ Step()

std::vector< std::shared_ptr< SamplingState > > GMHKernel::Step ( unsigned int const  t,
std::shared_ptr< SamplingState state 
)
overridevirtual

Sample GMHKernel::M states from the proposed states GMHKernel::proposedStates.

Parameters
[in]tThe current step in the MCMC chain
[in]stateThe current MCMC state
Returns
The accepted states

Reimplemented from muq::SamplingAlgorithms::MHKernel.

Definition at line 210 of file GMHKernel.cpp.

References muq::SamplingAlgorithms::TransitionKernel::comm, M, and SampleStationary().

Member Data Documentation

◆ M

const unsigned int muq::SamplingAlgorithms::GMHKernel::M
private

Number of accepted points (number of points added to the chain)

Definition at line 112 of file GMHKernel.h.

Referenced by SampleStationary(), and Step().

◆ N

const unsigned int muq::SamplingAlgorithms::GMHKernel::N
private

Number of proposals.

Definition at line 103 of file GMHKernel.h.

◆ Np1

const unsigned int muq::SamplingAlgorithms::GMHKernel::Np1
private

Number of proposals plus one.

Defined so we don't aways have to compute \(N+1\).

Definition at line 109 of file GMHKernel.h.

Referenced by AcceptanceDensity(), AcceptanceMatrix(), ComputeStationaryAcceptance(), SerialProposal(), and StationaryAcceptance().

◆ proposedStates

std::vector<std::shared_ptr<SamplingState> > muq::SamplingAlgorithms::GMHKernel::proposedStates
private

Proposed states.

Definition at line 118 of file GMHKernel.h.

Referenced by AcceptanceDensity(), SampleStationary(), and SerialProposal().

◆ stationaryAcceptance

Eigen::VectorXd muq::SamplingAlgorithms::GMHKernel::stationaryAcceptance
private

The cumulative stationary accepatnce probability.

Definition at line 115 of file GMHKernel.h.

Referenced by ComputeStationaryAcceptance(), SampleStationary(), and StationaryAcceptance().


The documentation for this class was generated from the following files: