MUQ  0.4.3
ModPieceWrapper.cpp
Go to the documentation of this file.
1 #include "AllClassWrappers.h"
2 
3 #include <pybind11/pybind11.h>
4 #include <pybind11/stl.h>
5 #include <pybind11/eigen.h>
6 #include "pybind11_json.h"
7 
11 #include "MUQ/Modeling/ModPiece.h"
18 #include "MUQ/Modeling/WorkGraph.h"
19 #include "MUQ/Modeling/SumPiece.h"
20 
21 #include <string>
22 
23 #include <functional>
24 #include <vector>
25 
26 using namespace muq::Modeling;
27 namespace py = pybind11;
28 
29 class PyModPieceTramp : public PyModPiece {
30 public:
31  /* Inherit the constructors */
33 
34  /* Trampoline (need one for each virtual function) */
35  void EvaluateImpl(std::vector<Eigen::VectorXd> const& input) override {
36  PYBIND11_OVERLOAD_PURE(void, PyModPiece, EvaluateImpl, input);
37  }
38 
39  void GradientImpl(unsigned int const outputDimWrt,
40  unsigned int const inputDimWrt,
41  std::vector<Eigen::VectorXd> const& input,
42  Eigen::VectorXd const& sensitivity) override {
43  PYBIND11_OVERLOAD(void, PyModPiece, GradientImpl, outputDimWrt, inputDimWrt, input, sensitivity);
44  }
45 
46  void JacobianImpl(unsigned int const outputDimWrt,
47  unsigned int const inputDimWrt,
48  std::vector<Eigen::VectorXd> const& input) override {
49  PYBIND11_OVERLOAD(void, PyModPiece, JacobianImpl, outputDimWrt, inputDimWrt, input);
50  }
51 
52  void ApplyJacobianImpl(unsigned int const outputDimWrt,
53  unsigned int const inputDimWrt,
54  std::vector<Eigen::VectorXd> const& input,
55  Eigen::VectorXd const& vec) override {
56  PYBIND11_OVERLOAD(void, PyModPiece, ApplyJacobianImpl, outputDimWrt, inputDimWrt, input, vec);
57  }
58 
59  void ApplyHessianImpl(unsigned int const outputDimWrt,
60  unsigned int const inputDimWrt1,
61  unsigned int const inputDimWrt2,
62  std::vector<Eigen::VectorXd> const& input,
63  Eigen::VectorXd const& sens,
64  Eigen::VectorXd const& vec) override {
65  PYBIND11_OVERLOAD(void, PyModPiece, ApplyHessianImpl, outputDimWrt, inputDimWrt1, inputDimWrt2, input, sens, vec);
66  }
67 };
68 
69 // The server itself does not do any Python calls; nevertheless, Python's global interpreter
70 // lock (GIL) remains active while entering the server loop, unless explicitly instructed not to.
71 // The server spawns threads for requests which may have to evaluate a Python model, in turn locking the GIL.
72 // We therefore have to unlock the GIL for the server loop, else we get a deadlock once
73 // a Python model is evaluated by the server.
74 void serveModPieceWithoutGIL(std::shared_ptr<ModPiece> modPiece, std::string name, std::string host, int port) {
75  Py_BEGIN_ALLOW_THREADS
76  muq::Modeling::serveModPiece(modPiece, name, host, port);
77  Py_END_ALLOW_THREADS
78 }
79 
80 class Publicist : public PyModPiece {
81 public:
82  // Expose protected functions
88 
89  // Expose protected member variables
90  using PyModPiece::outputs;
95 };
96 
98 {
99  // Define some functions from the WorkPiece base class
100  py::class_<ModPiece, WorkPiece, std::shared_ptr<ModPiece>> mp(m, "ModPiece");
101  mp
102  .def("Evaluate", (std::vector<Eigen::VectorXd> const& (ModPiece::*)(std::vector<Eigen::VectorXd> const&)) &ModPiece::Evaluate)
103  .def("Evaluate", (std::vector<Eigen::VectorXd> const& (ModPiece::*)()) &ModPiece::Evaluate)
104  .def_readonly("inputSizes", &ModPiece::inputSizes)
105  .def_readonly("outputSizes", &ModPiece::outputSizes)
106  .def("GetRunTime", &ModPiece::GetRunTime)
107  .def("ResetCallTime", &ModPiece::ResetCallTime)
108  .def("GetNumCalls", &ModPiece::GetNumCalls)
109  .def("Gradient", (Eigen::VectorXd const& (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &ModPiece::Gradient)
110  .def("Jacobian", (Eigen::MatrixXd const& (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&)) &ModPiece::Jacobian)
111  .def("ApplyJacobian", (Eigen::VectorXd const& (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &ModPiece::ApplyJacobian)
112  .def("GradientByFD", (Eigen::VectorXd (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &ModPiece::GradientByFD)
113  .def("JacobianByFD", (Eigen::MatrixXd (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&)) &ModPiece::JacobianByFD)
114  .def("ApplyJacobianByFD", (Eigen::VectorXd (ModPiece::*)(unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &ModPiece::ApplyJacobianByFD)
115  .def("ApplyHessian", (Eigen::VectorXd const& (ModPiece::*)(unsigned int, unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&, Eigen::VectorXd const&)) &ModPiece::ApplyHessian)
116  .def("ApplyHessianByFD", (Eigen::VectorXd (ModPiece::*)(unsigned int, unsigned int, unsigned int, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&, Eigen::VectorXd const&)) &ModPiece::ApplyHessianByFD)
117  .def("EnableCache", &ModPiece::EnableCache)
118  .def("DisableCache", &ModPiece::DisableCache)
119  .def("CacheStatus", &ModPiece::CacheStatus)
120  .def("SetWarnLevel", &ModPiece::SetWarnLevel);
121 
122  py::class_<PyModPiece, PyModPieceTramp, ModPiece, WorkPiece, std::shared_ptr<PyModPiece> > pymp(m, "PyModPiece");
123  pymp
124  .def(py::init<Eigen::VectorXi const&, Eigen::VectorXi const&>())
125  .def("EvaluateImpl", (void (PyModPiece::*)(std::vector<Eigen::VectorXd> const&)) &Publicist::EvaluateImpl)
126  .def_readwrite("outputs", &Publicist::outputs)
127  .def("GradientImpl", (void (PyModPiece::*)(unsigned int const, unsigned int const, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &Publicist::GradientImpl)
128  .def_readwrite("gradient", &Publicist::gradient)
129  .def("JacobianImpl", (void (PyModPiece::*)(unsigned int const, unsigned int const, std::vector<Eigen::VectorXd> const&)) &Publicist::JacobianImpl)
130  .def_readwrite("jacobian", &Publicist::jacobian)
131  .def("ApplyJacobianImpl", (void (PyModPiece::*)(unsigned int const, unsigned int const, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&)) &Publicist::ApplyJacobianImpl)
132  .def_readwrite("jacobianAction", &Publicist::jacobianAction)
133  .def("ApplyHessianImpl", (void (PyModPiece::*)(unsigned int const, unsigned int const, unsigned int const, std::vector<Eigen::VectorXd> const&, Eigen::VectorXd const&, Eigen::VectorXd const&)) &Publicist::ApplyHessianImpl)
134  .def_readwrite("hessAction", &Publicist::hessAction);
135 
136  py::class_<OneStepCachePiece, ModPiece, WorkPiece, std::shared_ptr<OneStepCachePiece>> ocp(m, "OneStepCachePiece");
137  ocp
138  .def(py::init<std::shared_ptr<ModPiece>>())
139  .def("HitRatio", &OneStepCachePiece::HitRatio);
140 
141  py::class_<UMBridgeModPiece, ModPiece, WorkPiece, std::shared_ptr<UMBridgeModPiece>> hmp(m, "UMBridgeModPiece");
142  hmp
143  .def(py::init( [](std::string host, std::string name) {return new UMBridgeModPiece(host, name); }))
144  .def(py::init( [](std::string host, std::string name, py::dict config) {return new UMBridgeModPiece(host, name, config); }));
145 
146  m.def("serveModPiece", &serveModPieceWithoutGIL);
147 
148  py::class_<ConstantVector, ModPiece, WorkPiece, std::shared_ptr<ConstantVector>> cv(m, "ConstantVector");
149  cv
150  .def(py::init<Eigen::VectorXd const&>())
151  .def("SetValue", &ConstantVector::SetValue);
152 
153  py::class_<ReplicateOperator, ModPiece, WorkPiece, std::shared_ptr<ReplicateOperator>> ro(m, "ReplicateOperator");
154  ro
155  .def(py::init<unsigned int, unsigned int>());
156 
157  py::class_<SumPiece, ModPiece, WorkPiece, std::shared_ptr<SumPiece>>(m, "SumPiece")
158  .def(py::init<unsigned int, unsigned int>(), py::arg("dim"),py::arg("numInputs")=2);
159 
160 
161  py::class_<ModGraphPiece, ModPiece, WorkPiece, std::shared_ptr<ModGraphPiece>> mgp(m, "ModGraphPiece");
162  mgp
163  .def("GetGraph", &ModGraphPiece::GetGraph)
164  .def("GetConstantPieces", &ModGraphPiece::GetConstantPieces)
165  .def("GradientGraph", &ModGraphPiece::GradientGraph)
166  .def("JacobianGraph", &ModGraphPiece::JacobianGraph);
167 
168  py::class_<MultiLogisticLikelihood, ModPiece, WorkPiece, std::shared_ptr<MultiLogisticLikelihood>> mll(m, "MultiLogisticLikelihood");
169  mll
170  .def(py::init<unsigned int, Eigen::VectorXi>());
171 
172  py::class_<SplitVector, ModPiece, WorkPiece, std::shared_ptr<SplitVector>>(m,"SplitVector")
173  .def(py::init<Eigen::VectorXi const&, Eigen::VectorXi const&, unsigned int const>())
174  .def("StartIndices", &SplitVector::StartIndices);
175 
176 }
void serveModPieceWithoutGIL(std::shared_ptr< ModPiece > modPiece, std::string name, std::string host, int port)
void SetValue(Eigen::VectorXd const &valIn)
Set the outputs.
std::vector< std::shared_ptr< ConstantVector > > GetConstantPieces()
Definition: ModGraphPiece.h:44
std::shared_ptr< WorkGraph > GetGraph()
Definition: ModGraphPiece.h:42
std::shared_ptr< ModGraphPiece > JacobianGraph(unsigned int const outputDimWrt, unsigned int const inputDimWrt)
std::shared_ptr< ModGraphPiece > GradientGraph(unsigned int const outputDimWrt, unsigned int const inputDimWrt)
Provides an abstract interface for defining vector-valued model components.
Definition: ModPiece.h:148
Eigen::MatrixXd jacobian
Definition: ModPiece.h:506
virtual Eigen::VectorXd GradientByFD(unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity)
Definition: ModPiece.cpp:111
virtual Eigen::VectorXd const & Gradient(unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity)
Compute the Gradient .
Definition: ModPiece.cpp:83
bool CacheStatus() const
Definition: ModPiece.h:461
virtual Eigen::VectorXd ApplyHessianByFD(unsigned int const outWrt, unsigned int const inWrt1, unsigned int const inWrt2, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec)
Definition: ModPiece.cpp:448
const Eigen::VectorXi inputSizes
Definition: ModPiece.h:469
virtual Eigen::MatrixXd const & Jacobian(unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input)
Compute the Jacobian of this ModPiece.
Definition: ModPiece.cpp:128
Eigen::MatrixXd ApplyJacobianByFD(unsigned int outWrt, unsigned int inWrt, Args const &... args)
Definition: ModPiece.h:303
std::vector< Eigen::VectorXd > outputs
Definition: ModPiece.h:503
Eigen::VectorXd gradient
Definition: ModPiece.h:504
virtual double GetRunTime(const std::string &method="Evaluate") const override
Get the average run time for one of the implemented methods.
Definition: ModPiece.cpp:497
virtual void SetWarnLevel(unsigned int newLevel)
Definition: ModPiece.h:469
virtual Eigen::VectorXd const & ApplyJacobian(unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec)
Apply the Jacobian of this ModPiece to a vector.
Definition: ModPiece.cpp:153
Eigen::VectorXd const & ApplyHessian(unsigned int outWrt, unsigned int inWrt1, unsigned int inWrt2, Eigen::VectorXd const &last, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec)
Definition: ModPiece.h:245
const Eigen::VectorXi outputSizes
Definition: ModPiece.h:472
Eigen::MatrixXd JacobianByFD(unsigned int outWrt, unsigned int inWrt, Args const &... args)
Definition: ModPiece.h:297
virtual unsigned long int GetNumCalls(const std::string &method="Evaluate") const override
get the number of times one of the implemented methods has been called.
Definition: ModPiece.cpp:535
virtual void ResetCallTime() override
Resets the number of call and times.
Definition: ModPiece.cpp:520
Eigen::VectorXd hessAction
Definition: ModPiece.h:507
Eigen::VectorXd jacobianAction
Definition: ModPiece.h:505
virtual void GradientImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) override
Definition: PyModPiece.cpp:13
virtual void ApplyJacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) override
Definition: PyModPiece.cpp:42
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &input) override
Definition: PyModPiece.cpp:9
virtual void ApplyHessianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt1, unsigned int const inputDimWrt2, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sens, Eigen::VectorXd const &vec)
Definition: PyModPiece.cpp:67
virtual void JacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input) override
Definition: PyModPiece.cpp:28
PyModPiece(Eigen::VectorXi const &inputSizes, Eigen::VectorXi const &outputSizes)
Definition: PyModPiece.cpp:5
Eigen::VectorXi StartIndices() const
Definition: SplitVector.h:64
A ModPiece connecting to a model via the UM-Bridge HTTP protocol.
std::vector< boost::any > const & Evaluate()
Evaluate this muq::Modeling::WorkPiece in the case that there are no inputs.
Definition: WorkPiece.cpp:222
void ModPieceWrapper(pybind11::module &m)
void serveModPiece(std::shared_ptr< ModPiece > modPiece, std::string name, std::string host, int port)
Serve a ModPiece via network using UM-Bridge.