1 #include "AllClassWrappers.h"
3 #include <pybind11/pybind11.h>
4 #include <pybind11/stl.h>
5 #include <pybind11/eigen.h>
35 void EvaluateImpl(std::vector<Eigen::VectorXd>
const& input)
override {
36 PYBIND11_OVERLOAD_PURE(
void,
PyModPiece, EvaluateImpl, input);
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);
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);
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);
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);
75 Py_BEGIN_ALLOW_THREADS
100 py::class_<ModPiece, WorkPiece, std::shared_ptr<ModPiece>> mp(m,
"ModPiece");
102 .def(
"Evaluate", (std::vector<Eigen::VectorXd>
const& (
ModPiece::*)(std::vector<Eigen::VectorXd>
const&)) &
ModPiece::Evaluate)
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)
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)
122 py::class_<PyModPiece, PyModPieceTramp, ModPiece, WorkPiece, std::shared_ptr<PyModPiece> > pymp(m,
"PyModPiece");
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);
136 py::class_<OneStepCachePiece, ModPiece, WorkPiece, std::shared_ptr<OneStepCachePiece>> ocp(m,
"OneStepCachePiece");
138 .def(py::init<std::shared_ptr<ModPiece>>())
141 py::class_<UMBridgeModPiece, ModPiece, WorkPiece, std::shared_ptr<UMBridgeModPiece>> hmp(m,
"UMBridgeModPiece");
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); }));
148 py::class_<ConstantVector, ModPiece, WorkPiece, std::shared_ptr<ConstantVector>> cv(m,
"ConstantVector");
150 .def(py::init<Eigen::VectorXd const&>())
153 py::class_<ReplicateOperator, ModPiece, WorkPiece, std::shared_ptr<ReplicateOperator>> ro(m,
"ReplicateOperator");
155 .def(py::init<unsigned int, unsigned int>());
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);
161 py::class_<ModGraphPiece, ModPiece, WorkPiece, std::shared_ptr<ModGraphPiece>> mgp(m,
"ModGraphPiece");
168 py::class_<MultiLogisticLikelihood, ModPiece, WorkPiece, std::shared_ptr<MultiLogisticLikelihood>> mll(m,
"MultiLogisticLikelihood");
170 .def(py::init<unsigned int, Eigen::VectorXi>());
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>())
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()
std::shared_ptr< WorkGraph > GetGraph()
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.
virtual Eigen::VectorXd GradientByFD(unsigned int const outputDimWrt, unsigned int const inputDimWrt, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity)
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 .
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)
const Eigen::VectorXi inputSizes
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.
Eigen::MatrixXd ApplyJacobianByFD(unsigned int outWrt, unsigned int inWrt, Args const &... args)
std::vector< Eigen::VectorXd > outputs
virtual double GetRunTime(const std::string &method="Evaluate") const override
Get the average run time for one of the implemented methods.
virtual void SetWarnLevel(unsigned int newLevel)
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.
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)
const Eigen::VectorXi outputSizes
Eigen::MatrixXd JacobianByFD(unsigned int outWrt, unsigned int inWrt, Args const &... args)
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.
virtual void ResetCallTime() override
Resets the number of call and times.
Eigen::VectorXd hessAction
Eigen::VectorXd jacobianAction
virtual void GradientImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) override
virtual void ApplyJacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) override
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &input) override
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)
virtual void JacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input) override
PyModPiece(Eigen::VectorXi const &inputSizes, Eigen::VectorXi const &outputSizes)
Eigen::VectorXi StartIndices() const
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.
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.