1 #include "AllClassWrappers.h"
10 #include <pybind11/pybind11.h>
11 #include <pybind11/stl.h>
12 #include <pybind11/eigen.h>
24 py::class_<TransitionKernel, std::shared_ptr<TransitionKernel>> transKern(m,
"TransitionKernel");
34 py::class_<MHKernel, TransitionKernel, std::shared_ptr<MHKernel>> mhKern(m,
"MHKernel");
36 .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem) {
return new MHKernel(
ConvertDictToPtree(d), problem);}))
37 .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem, std::shared_ptr<MCMCProposal> proposalIn) {
return new MHKernel(
ConvertDictToPtree(d), problem, proposalIn);}))
43 py::class_<DRKernel, TransitionKernel, std::shared_ptr<DRKernel>>(m,
"DRKernel")
44 .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem) {
return new DRKernel(
ConvertDictToPtree(d), problem);}))
45 .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem, std::vector<std::shared_ptr<MCMCProposal>> proposals, std::vector<double> scales) {
return new DRKernel(
ConvertDictToPtree(d), problem, proposals, scales);}))
46 .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem, std::vector<std::shared_ptr<MCMCProposal>> proposals) {
return new DRKernel(
ConvertDictToPtree(d), problem, proposals);}))
54 py::class_<DILIKernel, TransitionKernel, std::shared_ptr<DILIKernel>>(m,
"DILIKernel")
55 .def(py::init( [](py::dict d,
57 .def(py::init( [](py::dict d,
58 std::shared_ptr<AbstractSamplingProblem> problem,
59 std::shared_ptr<muq::Modeling::GaussianBase>
const& prior,
60 std::shared_ptr<muq::Modeling::ModPiece>
const& noiseMod,
61 std::shared_ptr<muq::Modeling::ModPiece>
const& likelihood){
return new DILIKernel(
ConvertDictToPtree(d),problem,prior,noiseMod, likelihood);}))
62 .def(py::init( [](py::dict d,
63 std::shared_ptr<AbstractSamplingProblem> problem,
64 std::shared_ptr<muq::Modeling::GaussianBase>
const& prior,
66 .def(py::init( [](py::dict d,
67 std::shared_ptr<AbstractSamplingProblem> problem,
68 Eigen::VectorXd
const& vals,
70 .def(py::init( [](py::dict d,
71 std::shared_ptr<AbstractSamplingProblem> problem,
72 std::shared_ptr<muq::Modeling::GaussianBase>
const& prior,
73 std::shared_ptr<muq::Modeling::ModPiece>
const& noiseMod,
74 std::shared_ptr<muq::Modeling::ModPiece>
const& likelihood,
75 Eigen::VectorXd
const& vals,
77 .def(py::init( [](py::dict d,
78 std::shared_ptr<AbstractSamplingProblem> problem,
79 std::shared_ptr<muq::Modeling::GaussianBase>
const& prior,
80 std::shared_ptr<muq::Modeling::ModPiece>
const& likelihood,
81 Eigen::VectorXd
const& vals,
An implementation of the Dimension Independent Likelihood Informed subspace (DILI) MCMC sampler.
virtual std::shared_ptr< TransitionKernel > LISKernel()
static std::shared_ptr< muq::Modeling::ModPiece > CreateLikelihood(std::shared_ptr< muq::Modeling::ModPiece > const &forwardModel, std::shared_ptr< muq::Modeling::ModPiece > const &noiseDensity)
Eigen::VectorXd FromLIS(Eigen::VectorXd const &r) const
Eigen::VectorXd const & LISVals() const
static std::shared_ptr< muq::Modeling::ModPiece > ExtractLikelihood(std::shared_ptr< AbstractSamplingProblem > const &problem, std::string const &nodeName)
unsigned int LISDim() const
Eigen::MatrixXd const & LISVecs() const
static std::shared_ptr< muq::Modeling::ModPiece > ExtractForwardModel(std::shared_ptr< muq::Modeling::ModPiece > const &likelihoodIn)
virtual std::shared_ptr< TransitionKernel > CSKernel()
static std::shared_ptr< muq::Modeling::ModPiece > ExtractNoiseModel(std::shared_ptr< muq::Modeling::ModPiece > const &likelihood)
Eigen::VectorXd ToLIS(Eigen::VectorXd const &x) const
void CreateLIS(std::vector< Eigen::VectorXd > const &currState)
Eigen::VectorXd ToCS(Eigen::VectorXd const &x) const
An implementation of the delayed rejection kernel.
double EvaluateProposal(unsigned int stage, std::shared_ptr< SamplingState > const &x, std::shared_ptr< SamplingState > const &y) const
virtual Eigen::VectorXd AcceptanceRates() const
virtual std::vector< std::shared_ptr< MCMCProposal > > Proposals()
virtual std::vector< std::shared_ptr< SamplingState > > Step(unsigned int const t, std::shared_ptr< SamplingState > prevState) override
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.
std::shared_ptr< SamplingState > SampleProposal(unsigned int stage, std::shared_ptr< SamplingState > const &state) const
An implementation of the standard Metropolis-Hastings transition kernel.
virtual double AcceptanceRate() const
virtual std::vector< std::shared_ptr< SamplingState > > Step(unsigned int const t, std::shared_ptr< SamplingState > prevState) override
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.
virtual void PreStep(unsigned int const t, std::shared_ptr< SamplingState > state)
Allow the kernel to preprocess the current step.
virtual int GetBlockInd() const
static std::shared_ptr< TransitionKernel > Construct(boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > problem)
Static constructor for the transition kernel.
virtual std::vector< std::shared_ptr< SamplingState > > Step(unsigned int const t, std::shared_ptr< SamplingState > prevState)=0
virtual void SetBlockInd(int newBlockInd)
virtual void PostStep(unsigned int const t, std::vector< std::shared_ptr< SamplingState >> const &state)
Allow the kernel to adapt given a new state.
void KernelWrapper(pybind11::module &m)
boost::property_tree::ptree ConvertDictToPtree(pybind11::dict dict)