A kernel for the generalized Metropolis-Hastings kernel. More...
#include <GMHKernel.h>
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< MCMCProposal > | Proposal () |
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, TransitionKernelConstructor > | TransitionKernelMap |
Static Public Member Functions inherited from muq::SamplingAlgorithms::TransitionKernel | |
static std::shared_ptr< TransitionKernel > | Construct (boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem) |
Static constructor for the transition kernel. More... | |
static std::shared_ptr< TransitionKernelMap > | GetTransitionKernelMap () |
muq::SamplingAlgorithms::GMHKernel::GMHKernel | ( | boost::property_tree::ptree const & | pt, |
std::shared_ptr< AbstractSamplingProblem > | problem | ||
) |
[in] | pt | Options for this kenel and the standard Metropolis-Hastings kernel |
[in] | problem | The problem we want to sample |
muq::SamplingAlgorithms::GMHKernel::GMHKernel | ( | boost::property_tree::ptree const & | pt, |
std::shared_ptr< AbstractSamplingProblem > | problem, | ||
std::shared_ptr< MCMCProposal > | proposalIn | ||
) |
[in] | pt | Options for this kenel and the standard Metropolis-Hastings kernel |
[in] | problem | The problem we want to sample |
[in] | proposalIn | The proposal for the MCMC chain |
|
virtual |
Definition at line 56 of file GMHKernel.cpp.
|
private |
Compute the stationary transition density.
[in] | R | The 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().
|
private |
Definition at line 154 of file GMHKernel.cpp.
References Np1.
Referenced by ComputeStationaryAcceptance().
|
private |
Compute the cumulative acceptance density.
[in] | R | The log-target plus the log-proposals |
Definition at line 168 of file GMHKernel.cpp.
References AcceptanceMatrix(), Np1, and stationaryAcceptance.
Referenced by AcceptanceDensity().
std::shared_ptr< parcer::Communicator > GMHKernel::GetCommunicator | ( | ) | const |
Definition at line 228 of file GMHKernel.cpp.
References muq::SamplingAlgorithms::TransitionKernel::comm.
|
private |
Propose \(N\) points in parallel and evaluate the log target.
[in] | t | The current step in the MCMC chain |
[in] | state | The current point |
Definition at line 84 of file GMHKernel.cpp.
References muq::SamplingAlgorithms::TransitionKernel::comm, and SerialProposal().
Referenced by PreStep().
|
overridevirtual |
Propose GMHKernel::N points and compute the cumulative distribution of the stationary distribution for the acceptance probability (GMHKernel::proposedStates and GMHKernel::stationaryAcceptance, respectively)
[in] | t | The current step in the MCMC chain |
[in] | state | The current MCMC state |
Reimplemented from muq::SamplingAlgorithms::TransitionKernel.
Definition at line 183 of file GMHKernel.cpp.
References muq::SamplingAlgorithms::TransitionKernel::comm, ParallelProposal(), and SerialProposal().
std::vector< std::shared_ptr< SamplingState > > GMHKernel::SampleStationary | ( | ) | const |
Sample the stationary distribution.
Definition at line 197 of file GMHKernel.cpp.
References M, proposedStates, and stationaryAcceptance.
Referenced by Step().
|
private |
Propose \(N\) points in serial and evaluate the log target.
[in] | t | The current step in the MCMC chain |
[in] | state | The 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().
Eigen::VectorXd GMHKernel::StationaryAcceptance | ( | ) | const |
Get the cumulative stationary acceptance probability.
Must be called after GMHKernel::PreStep.
Definition at line 220 of file GMHKernel.cpp.
References Np1, and stationaryAcceptance.
|
overridevirtual |
Sample GMHKernel::M states from the proposed states GMHKernel::proposedStates.
[in] | t | The current step in the MCMC chain |
[in] | state | The current MCMC state |
Reimplemented from muq::SamplingAlgorithms::MHKernel.
Definition at line 210 of file GMHKernel.cpp.
References muq::SamplingAlgorithms::TransitionKernel::comm, M, and SampleStationary().
|
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().
|
private |
Number of proposals.
Definition at line 103 of file GMHKernel.h.
|
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().
|
private |
Proposed states.
Definition at line 118 of file GMHKernel.h.
Referenced by AcceptanceDensity(), SampleStationary(), and SerialProposal().
|
private |
The cumulative stationary accepatnce probability.
Definition at line 115 of file GMHKernel.h.
Referenced by ComputeStationaryAcceptance(), SampleStationary(), and StationaryAcceptance().