MUQ  0.4.3
MCMCWrapper.cpp
Go to the documentation of this file.
1 #include "AllClassWrappers.h"
2 
3 #include "MUQ/config.h"
4 
11 
14 
15 #include <pybind11/pybind11.h>
16 #include <pybind11/stl.h>
17 #include <pybind11/eigen.h>
18 #include <pybind11/iostream.h>
19 #include <pybind11/numpy.h>
20 
21 #include <string>
22 
23 #include <functional>
24 #include <vector>
25 
26 using namespace muq::SamplingAlgorithms;
27 using namespace muq::Utilities;
28 namespace py = pybind11;
29 
34 
35 using namespace muq::Modeling;
36 
37 
38 void PythonBindings::MCMCWrapper(py::module &m) {
39 
40  py::class_<SingleChainMCMC, std::shared_ptr<SingleChainMCMC>> singleMCMC(m, "SingleChainMCMC");
41  singleMCMC
42  .def(py::init( [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem) {return new SingleChainMCMC(ConvertDictToPtree(d), problem);}))
43  .def(py::init( [](py::dict d, std::vector<std::shared_ptr<TransitionKernel>> kernels) {return new SingleChainMCMC(ConvertDictToPtree(d), kernels);}))
44  .def("SetState", (void (SingleChainMCMC::*)(std::shared_ptr<SamplingState> const&)) &SingleChainMCMC::SetState)
45  .def("SetState", (void (SingleChainMCMC::*)(std::vector<Eigen::VectorXd> const&)) &SingleChainMCMC::SetState)
46  .def("Kernels", &SingleChainMCMC::Kernels)
47  .def("Run", (std::shared_ptr<MarkovChain> (SingleChainMCMC::*)(std::vector<Eigen::VectorXd> const&)) &SingleChainMCMC::Run, py::call_guard<py::scoped_ostream_redirect, py::scoped_estream_redirect>())
48  .def("AddNumSamps", &SingleChainMCMC::AddNumSamps)
49  .def("NumSamps", &SingleChainMCMC::NumSamps)
50  .def("TotalTime", &SingleChainMCMC::TotalTime)
51  .def("GetSamples", &SingleChainMCMC::GetSamples)
52  .def("GetQOIs", &SingleChainMCMC::GetQOIs);
53 
54  py::class_<ParallelTempering, std::shared_ptr<ParallelTempering>>(m,"ParallelTempering")
55  .def(py::init( [](py::dict d, std::shared_ptr<InferenceProblem> problem) {return new ParallelTempering(ConvertDictToPtree(d), problem);}))
56  .def(py::init( [](py::dict d, Eigen::VectorXd const& invTemps, std::vector<std::shared_ptr<TransitionKernel>> kerns) {return new ParallelTempering(ConvertDictToPtree(d), invTemps, kerns);}))
57  .def(py::init( [](py::dict d, Eigen::VectorXd const& invTemps, std::vector<std::vector<std::shared_ptr<TransitionKernel>>> kerns) {return new ParallelTempering(ConvertDictToPtree(d), invTemps, kerns);}))
58  .def("SetState", (void (ParallelTempering::*)(std::vector<std::shared_ptr<SamplingState>> const&)) &ParallelTempering::SetState)
59  .def("SetState", (void (ParallelTempering::*)(std::vector<Eigen::VectorXd> const&)) &ParallelTempering::SetState)
60  .def("SetState", (void (ParallelTempering::*)(std::vector<std::vector<Eigen::VectorXd>> const&)) &ParallelTempering::SetState)
61  .def("GetInverseTemp", &ParallelTempering::GetInverseTemp)
62  .def("Kernels", &ParallelTempering::Kernels)
63  .def("Run", (std::shared_ptr<MarkovChain> (ParallelTempering::*)(Eigen::VectorXd const&)) &ParallelTempering::Run)
64  .def("Run", (std::shared_ptr<MarkovChain> (ParallelTempering::*)(std::vector<Eigen::VectorXd> const&)) &ParallelTempering::Run)
65  .def("Run", (std::shared_ptr<MarkovChain> (ParallelTempering::*)(std::vector<std::vector<Eigen::VectorXd>> const&)) &ParallelTempering::Run)
66  .def("AddNumSamps", &ParallelTempering::AddNumSamps)
67  .def("NumSamps", &ParallelTempering::NumSamps)
68  .def("GetSamples", &ParallelTempering::GetSamples)
69  .def("GetQOISs", &ParallelTempering::GetQOIs)
70  .def_readonly("numTemps", &ParallelTempering::numTemps);
71 
72 
73  py::class_<MIMCMCBox, std::shared_ptr<MIMCMCBox>> multiindexMCMCBox(m, "MIMCMCBox");
74  multiindexMCMCBox
75  .def("FinestChain", &MIMCMCBox::FinestChain)
76  .def("GetChain", &MIMCMCBox::GetChain)
77  .def("GetBoxIndices", &MIMCMCBox::GetBoxIndices)
78  .def("GetHighestIndex", &MIMCMCBox::GetHighestIndex);
79 
80 
81  py::class_<MIMCMC, std::shared_ptr<MIMCMC>> multiindexMCMC(m, "MIMCMC");
82  multiindexMCMC
83  .def(py::init( [](py::dict d, Eigen::VectorXd startingPoint, std::vector<std::shared_ptr<AbstractSamplingProblem>> const& problems) {return new MIMCMC(ConvertDictToPtree(d), startingPoint, problems); }))
84  .def(py::init( [](py::dict d, Eigen::VectorXd startingPoint, std::vector<std::shared_ptr<ModPiece>> const& models) {return new MIMCMC(ConvertDictToPtree(d), startingPoint, models); }))
85  .def(py::init( [](py::dict d, Eigen::VectorXd startingPoint, std::vector<std::shared_ptr<AbstractSamplingProblem>> const& problems, std::shared_ptr<MultiIndexSet> const& indices) {return new MIMCMC(ConvertDictToPtree(d), startingPoint, problems, indices); }))
86  .def(py::init( [](py::dict d, Eigen::VectorXd startingPoint, std::vector<std::shared_ptr<ModPiece>> const& models, std::shared_ptr<MultiIndexSet> const& indices) {return new MIMCMC(ConvertDictToPtree(d), startingPoint, models, indices); }))
87  .def("Run", &MIMCMC::Run, py::call_guard<py::scoped_ostream_redirect, py::scoped_estream_redirect>())
88  .def("GetSamples", &MIMCMC::GetSamples)
89  .def("GetQOIs", &MIMCMC::GetQOIs)
90  .def("GetIndices", &MIMCMC::GetIndices)
91  .def("GetMIMCMCBox", &MIMCMC::GetMIMCMCBox);
92 
93  py::class_<MCMCFactory, std::shared_ptr<MCMCFactory>> fact(m, "MCMCFactory");
94  fact
95  .def_static("CreateSingleChain", [](py::dict d, std::shared_ptr<AbstractSamplingProblem> problem) {return MCMCFactory::CreateSingleChain(ConvertDictToPtree(d), problem);},
96  py::call_guard<py::scoped_ostream_redirect,py::scoped_estream_redirect>() );
97 
98 
99 }
static std::shared_ptr< SingleChainMCMC > CreateSingleChain(boost::property_tree::ptree pt, std::shared_ptr< AbstractSamplingProblem > problem)
Definition: MCMCFactory.cpp:8
std::shared_ptr< SingleChainMCMC > FinestChain()
Definition: MIMCMCBox.cpp:282
std::shared_ptr< MultiIndexSet > GetBoxIndices()
Definition: MIMCMCBox.cpp:287
std::shared_ptr< MultiIndex > GetHighestIndex()
Definition: MIMCMCBox.cpp:107
std::shared_ptr< SingleChainMCMC > GetChain(std::shared_ptr< MultiIndex > index)
Definition: MIMCMCBox.cpp:291
Multiindex MCMC method.
Definition: MIMCMC.h:21
std::shared_ptr< MultiIndexSet > GetIndices()
Get set of indices of boxes set up by the method.
Definition: MIMCMC.cpp:79
virtual std::shared_ptr< MultiIndexEstimator > GetQOIs() const
Definition: MIMCMC.cpp:48
std::shared_ptr< MIMCMCBox > GetMIMCMCBox(std::shared_ptr< MultiIndex > index)
Definition: MIMCMC.cpp:64
virtual std::shared_ptr< MultiIndexEstimator > Run()
Definition: MIMCMC.cpp:52
virtual std::shared_ptr< MultiIndexEstimator > GetSamples() const
Definition: MIMCMC.cpp:45
Defines an MCMC sampler with multiple chains running on problems with different temperatues.
const unsigned int numTemps
Number of temperatures in the temperature schedule.
std::shared_ptr< MarkovChain > GetQOIs() const
double GetInverseTemp(unsigned int chainInd) const
void SetState(std::vector< std::shared_ptr< SamplingState >> const &x0)
Set the state of the MCMC chain.
void AddNumSamps(unsigned int numNewSamps)
std::vector< std::shared_ptr< TransitionKernel > > const & Kernels(unsigned int chainInd) const
std::shared_ptr< MarkovChain > Run()
std::shared_ptr< MarkovChain > GetSamples() const
Defines an MCMC sampler with a single chain.
virtual std::shared_ptr< MarkovChain > GetSamples() const
virtual unsigned int NumSamps() const
std::shared_ptr< MarkovChain > Run(Args const &... args)
virtual void SetState(std::shared_ptr< SamplingState > const &x0)
Set the state of the MCMC chain.
virtual std::vector< std::shared_ptr< TransitionKernel > > & Kernels()
virtual void AddNumSamps(unsigned int numNewSamps)
std::shared_ptr< MarkovChain > GetQOIs() const
void MCMCWrapper(pybind11::module &m)
boost::property_tree::ptree ConvertDictToPtree(pybind11::dict dict)