MUQ  0.4.3
KernelWrapper.cpp
Go to the documentation of this file.
1 #include "AllClassWrappers.h"
2 
7 
9 
10 #include <pybind11/pybind11.h>
11 #include <pybind11/stl.h>
12 #include <pybind11/eigen.h>
13 
14 #include <string>
15 
16 #include <functional>
17 #include <vector>
18 
19 using namespace muq::SamplingAlgorithms;
20 using namespace muq::Utilities;
21 namespace py = pybind11;
22 
23 void PythonBindings::KernelWrapper(py::module &m) {
24  py::class_<TransitionKernel, std::shared_ptr<TransitionKernel>> transKern(m, "TransitionKernel");
25  transKern
26  .def_static("Construct", [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem)->std::shared_ptr<TransitionKernel>{return TransitionKernel::Construct(ConvertDictToPtree(d), problem);})
27  .def("PreStep", &TransitionKernel::PreStep)
28  .def("PostStep", &TransitionKernel::PostStep)
29  .def("Step", &TransitionKernel::Step)
30  .def("GetBlockInd", &TransitionKernel::GetBlockInd)
31  .def("SetBlockInd", &TransitionKernel::SetBlockInd);
32 
33 
34  py::class_<MHKernel, TransitionKernel, std::shared_ptr<MHKernel>> mhKern(m, "MHKernel");
35  mhKern
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);}))
38  .def("Proposal", &MHKernel::Proposal)
39  .def("PostStep", &MHKernel::PostStep)
40  .def("Step", &MHKernel::Step)
41  .def("AcceptanceRate", &MHKernel::AcceptanceRate);
42 
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);}))
47  .def("PostStep", &DRKernel::PostStep)
48  .def("Step", &DRKernel::Step)
49  .def("Proposals",&DRKernel::Proposals)
50  .def("AcceptanceRates",&DRKernel::AcceptanceRates)
51  .def("SampleProposal", &DRKernel::SampleProposal)
52  .def("EvaluateProposal", &DRKernel::EvaluateProposal);
53 
54  py::class_<DILIKernel, TransitionKernel, std::shared_ptr<DILIKernel>>(m,"DILIKernel")
55  .def(py::init( [](py::dict d,
56  std::shared_ptr<AbstractSamplingProblem> problem){return new DILIKernel(ConvertDictToPtree(d),problem);}))
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,
65  std::shared_ptr<muq::Modeling::ModPiece> const& likelihood){return new DILIKernel(ConvertDictToPtree(d),problem,prior,likelihood);}))
66  .def(py::init( [](py::dict d,
67  std::shared_ptr<AbstractSamplingProblem> problem,
68  Eigen::VectorXd const& vals,
69  Eigen::MatrixXd const& vecs){return new DILIKernel(ConvertDictToPtree(d),problem, vals, vecs);}))
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,
76  Eigen::MatrixXd const& vecs){return new DILIKernel(ConvertDictToPtree(d),problem,prior,noiseMod, likelihood, vals, vecs);}))
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,
82  Eigen::MatrixXd const& vecs){return new DILIKernel(ConvertDictToPtree(d),problem,prior,likelihood,vals,vecs);}))
83  .def("LISKernel", &DILIKernel::LISKernel)
84  .def("CSKernel", &DILIKernel::CSKernel)
85  .def_static("ExtractLikelihood", &DILIKernel::ExtractLikelihood)
86  .def_static("ExtractNoiseModel", &DILIKernel::ExtractNoiseModel)
87  .def_static("ExtractForwardModel", &DILIKernel::ExtractForwardModel)
88  .def_static("CreateLikelihood", &DILIKernel::CreateLikelihood)
89  .def("LISVecs", &DILIKernel::LISVecs)
90  .def("LISVals", &DILIKernel::LISVals)
91  .def("LISDim", &DILIKernel::LISDim)
92  .def("CreateLIS", &DILIKernel::CreateLIS)
93  .def("ToLIS", &DILIKernel::ToLIS)
94  .def("FromLIS", &DILIKernel::FromLIS)
95  .def("ToCS", &DILIKernel::ToCS);
96 }
An implementation of the Dimension Independent Likelihood Informed subspace (DILI) MCMC sampler.
Definition: DILIKernel.h:141
virtual std::shared_ptr< TransitionKernel > LISKernel()
Definition: DILIKernel.h:181
static std::shared_ptr< muq::Modeling::ModPiece > CreateLikelihood(std::shared_ptr< muq::Modeling::ModPiece > const &forwardModel, std::shared_ptr< muq::Modeling::ModPiece > const &noiseDensity)
Definition: DILIKernel.cpp:519
Eigen::VectorXd FromLIS(Eigen::VectorXd const &r) const
Definition: DILIKernel.cpp:223
Eigen::VectorXd const & LISVals() const
Definition: DILIKernel.h:214
static std::shared_ptr< muq::Modeling::ModPiece > ExtractLikelihood(std::shared_ptr< AbstractSamplingProblem > const &problem, std::string const &nodeName)
Definition: DILIKernel.cpp:498
Eigen::MatrixXd const & LISVecs() const
Definition: DILIKernel.h:213
static std::shared_ptr< muq::Modeling::ModPiece > ExtractForwardModel(std::shared_ptr< muq::Modeling::ModPiece > const &likelihoodIn)
Definition: DILIKernel.cpp:481
virtual std::shared_ptr< TransitionKernel > CSKernel()
Definition: DILIKernel.h:182
static std::shared_ptr< muq::Modeling::ModPiece > ExtractNoiseModel(std::shared_ptr< muq::Modeling::ModPiece > const &likelihood)
Definition: DILIKernel.cpp:471
Eigen::VectorXd ToLIS(Eigen::VectorXd const &x) const
Definition: DILIKernel.cpp:215
void CreateLIS(std::vector< Eigen::VectorXd > const &currState)
Definition: DILIKernel.cpp:372
Eigen::VectorXd ToCS(Eigen::VectorXd const &x) const
Definition: DILIKernel.cpp:229
An implementation of the delayed rejection kernel.
Definition: DRKernel.h:63
double EvaluateProposal(unsigned int stage, std::shared_ptr< SamplingState > const &x, std::shared_ptr< SamplingState > const &y) const
Definition: DRKernel.cpp:81
virtual Eigen::VectorXd AcceptanceRates() const
Definition: DRKernel.h:105
virtual std::vector< std::shared_ptr< MCMCProposal > > Proposals()
Definition: DRKernel.h:83
virtual std::vector< std::shared_ptr< SamplingState > > Step(unsigned int const t, std::shared_ptr< SamplingState > prevState) override
Definition: DRKernel.cpp:99
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.
Definition: DRKernel.cpp:63
std::shared_ptr< SamplingState > SampleProposal(unsigned int stage, std::shared_ptr< SamplingState > const &state) const
Definition: DRKernel.cpp:71
An implementation of the standard Metropolis-Hastings transition kernel.
Definition: MHKernel.h:22
virtual double AcceptanceRate() const
Definition: MHKernel.h:43
virtual std::vector< std::shared_ptr< SamplingState > > Step(unsigned int const t, std::shared_ptr< SamplingState > prevState) override
Definition: MHKernel.cpp:40
virtual std::shared_ptr< MCMCProposal > Proposal()
Definition: MHKernel.h:34
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.
Definition: MHKernel.cpp:36
virtual void PreStep(unsigned int const t, std::shared_ptr< SamplingState > state)
Allow the kernel to preprocess the current step.
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)