MUQ  0.4.3
PolynomialsWrapper.cpp
Go to the documentation of this file.
1 #include "AllClassWrappers.h"
2 
14 
15 #include <pybind11/pybind11.h>
16 #include <pybind11/stl.h>
17 #include <pybind11/eigen.h>
18 
19 #include <string>
20 
21 #include <functional>
22 #include <vector>
23 
25 using namespace muq::Modeling;
26 
27 namespace py = pybind11;
28 
30 {
31  py::class_<BasisExpansion, ModPiece, std::shared_ptr<BasisExpansion>> baseExp(m, "BasisExpansion");
32  baseExp
33  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&>())
34  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&, bool>())
35  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&, std::shared_ptr<muq::Utilities::MultiIndexSet>>())
36  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&, std::shared_ptr<muq::Utilities::MultiIndexSet>, bool>())
37  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&, std::shared_ptr<muq::Utilities::MultiIndexSet>, Eigen::MatrixXd const&>())
38  .def(py::init<std::vector<std::shared_ptr<IndexedScalarBasis>> const&, std::shared_ptr<muq::Utilities::MultiIndexSet>, Eigen::MatrixXd const&, bool>())
39  .def("NumTerms", &BasisExpansion::NumTerms)
40  .def("SecondDerivative", (Eigen::MatrixXd (BasisExpansion::*)(unsigned, unsigned, unsigned, Eigen::VectorXd const&, Eigen::MatrixXd const&)) &BasisExpansion::SecondDerivative)
41  .def("SecondDerivative", (Eigen::MatrixXd (BasisExpansion::*)(unsigned, unsigned, unsigned, Eigen::VectorXd const&)) &BasisExpansion::SecondDerivative)
42  .def("GetCoeffs", &BasisExpansion::GetCoeffs)
43  .def("SetCoeffs", &BasisExpansion::SetCoeffs)
44  .def("Multis", &BasisExpansion::Multis)
45  .def("BuildVandermonde", &BasisExpansion::BuildVandermonde)
46  .def("BuildDerivMatrix", &BasisExpansion::BuildDerivMatrix)
47  .def("ToHDF5", (void (BasisExpansion::*)(std::string, std::string) const) &BasisExpansion::ToHDF5, py::arg("filename"),py::arg("groupName")="/")
48  .def_static("FromHDF5", [](std::string filename, std::string dsetName){return BasisExpansion::FromHDF5(filename, dsetName);})
49  .def_static("FromHDF5", [](std::string filename){return BasisExpansion::FromHDF5(filename, "/");});
50 
51  py::class_<IndexedScalarBasis, std::shared_ptr<IndexedScalarBasis>> iSB(m, "IndexedScalarBasis");
52  iSB
53  .def_static("Construct", &IndexedScalarBasis::Construct)
54  .def("EvaluateAllTerms", &HermiteFunction::EvaluateAllTerms);
55 
56  py::class_<OrthogonalPolynomial, IndexedScalarBasis, std::shared_ptr<OrthogonalPolynomial>> orthPoly(m, "OrthogonalPolynomial");
57  orthPoly
58  .def_static("Construct", &OrthogonalPolynomial::Construct)
59  .def("BasisEvaluate", &OrthogonalPolynomial::BasisEvaluate)
60  .def("EvaluateAllTerms", &OrthogonalPolynomial::DerivativeEvaluate)
61  .def("Normalization", &OrthogonalPolynomial::Normalization);
62 
63  py::class_<HermiteFunction, IndexedScalarBasis, std::shared_ptr<HermiteFunction>> hermFunc(m, "HermiteFunction");
64  hermFunc
65  .def(py::init<>())
66  .def("BasisEvaluate", &HermiteFunction::BasisEvaluate)
67  .def("DerivativeEvaluate", &HermiteFunction::DerivativeEvaluate);
68 
69  py::class_<Jacobi, OrthogonalPolynomial, std::shared_ptr<Jacobi>> jacob(m, "Jacobi");
70  jacob
71  .def(py::init<>())
72  .def(py::init<const double>())
73  .def(py::init<const double, const double>())
74  .def("DerivativeEvaluate", &Jacobi::DerivativeEvaluate)
75  .def("Normalization", &Jacobi::Normalization);
76 
77  py::class_<Laguerre, OrthogonalPolynomial, std::shared_ptr<Laguerre>> laguerre(m, "Laguerre");
78  laguerre
79  .def(py::init<>())
80  .def(py::init<const double>())
81  .def("DerivativeEvaluate", &Laguerre::DerivativeEvaluate)
82  .def("Normalization", &Laguerre::Normalization);
83 
84  py::class_<Legendre, OrthogonalPolynomial, std::shared_ptr<Legendre>> legend(m, "Legendre");
85  legend
86  .def(py::init<>())
87  .def("DerivativeEvaluate", &Legendre::DerivativeEvaluate)
88  .def("Normalization", &Legendre::Normalization);
89 
90  py::class_<Monomial, IndexedScalarBasis, std::shared_ptr<Monomial>> mono(m, "Monomial");
91  mono
92  .def(py::init<>())
93  .def("BasisEvaluate", &Monomial::BasisEvaluate)
94  .def("DerivativeEvaluate", &Monomial::DerivativeEvaluate);
95 
96  py::class_<MonotoneExpansion, ModPiece, std::shared_ptr<MonotoneExpansion>> monotone(m, "MonotoneExpansion");
97  monotone
98  .def(py::init<std::shared_ptr<BasisExpansion>>())
99  .def(py::init<std::shared_ptr<BasisExpansion>, bool>())
100  .def(py::init<std::vector<std::shared_ptr<BasisExpansion>> const&, std::vector<std::shared_ptr<BasisExpansion>> const&>())
101  .def(py::init<std::vector<std::shared_ptr<BasisExpansion>> const&, std::vector<std::shared_ptr<BasisExpansion>> const&, bool>())
102  .def("NumTerms", &MonotoneExpansion::NumTerms)
103  .def("Head", &MonotoneExpansion::Head)
104  .def("EvaluateInverse", (Eigen::VectorXd (MonotoneExpansion::*)(Eigen::VectorXd const&) const) &MonotoneExpansion::EvaluateInverse)
105  .def("EvaluateInverse", (Eigen::VectorXd (MonotoneExpansion::*)(Eigen::VectorXd const&, Eigen::VectorXd const&) const) &MonotoneExpansion::EvaluateInverse)
106  .def("EvaluateForward", &MonotoneExpansion::EvaluateForward)
107  .def("GetCoeffs", &MonotoneExpansion::GetCoeffs)
108  .def("SetCoeffs", &MonotoneExpansion::SetCoeffs)
109  .def("LogDeterminant", (double (MonotoneExpansion::*)(Eigen::VectorXd const&)) &MonotoneExpansion::LogDeterminant)
110  .def("LogDeterminant", (double (MonotoneExpansion::*)(Eigen::VectorXd const&, Eigen::VectorXd const&)) &MonotoneExpansion::LogDeterminant)
111  .def("GradLogDeterminant", (Eigen::VectorXd (MonotoneExpansion::*)(Eigen::VectorXd const&)) &MonotoneExpansion::GradLogDeterminant)
112  .def("GradLogDeterminant", (Eigen::VectorXd (MonotoneExpansion::*)(Eigen::VectorXd const&, Eigen::VectorXd const&)) &MonotoneExpansion::GradLogDeterminant);
113 
114  py::class_<PhysicistHermite, OrthogonalPolynomial, std::shared_ptr<PhysicistHermite>> physHerm(m, "PhysicistHermite");
115  physHerm
116  .def(py::init<>())
117  .def("DerivativeEvaluate", &PhysicistHermite::DerivativeEvaluate)
118  .def("Normalization", &PhysicistHermite::Normalization);
119 
120  py::class_<ProbabilistHermite, OrthogonalPolynomial, std::shared_ptr<ProbabilistHermite>> probHerm(m, "ProbabilistHermite");
121  probHerm
122  .def(py::init<>())
123  .def("DerivativeEvaluate", &ProbabilistHermite::DerivativeEvaluate)
124  .def("Normalization", &ProbabilistHermite::Normalization);
125 }
Eigen::MatrixXd BuildDerivMatrix(Eigen::MatrixXd const &evalPts, int wrtDim) const
static std::shared_ptr< BasisExpansion > FromHDF5(std::string filename, std::string groupName="/")
Loads an expansion from an HDF5 file.
Eigen::MatrixXd BuildVandermonde(Eigen::MatrixXd const &evalPts) const
virtual unsigned NumTerms() const
const std::shared_ptr< muq::Utilities::MultiIndexSet > Multis() const
Eigen::MatrixXd GetCoeffs() const
Eigen::MatrixXd SecondDerivative(unsigned outputDim, unsigned wrtDim1, unsigned wrtDim2, Eigen::VectorXd const &evalPt, Eigen::MatrixXd const &coeffs)
void SetCoeffs(Eigen::MatrixXd const &allCoeffs)
virtual void ToHDF5(std::string filename, std::string groupName="/") const
Saves the expansion to group in an HDF5 file.
virtual double BasisEvaluate(int const order, double const x) const override
Evaluate the hermite function.
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
static std::shared_ptr< IndexedScalarBasis > Construct(std::string const &polyName)
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const =0
virtual Eigen::VectorXd EvaluateAllTerms(int const maxOrder, double const x) const
Evaluates all basis functions with order <= maxOrder.
virtual double Normalization(unsigned int polyOrder) const override
Definition: Jacobi.cpp:40
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
Definition: Jacobi.cpp:6
virtual double Normalization(unsigned int polyOrder) const override
Definition: Laguerre.cpp:34
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
Definition: Laguerre.cpp:6
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
Definition: Legendre.cpp:9
virtual double Normalization(unsigned int polyOrder) const override
Definition: Legendre.cpp:42
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
Definition: Monomial.cpp:9
virtual double BasisEvaluate(int const order, double const x) const override
Evaluate the specific basis type (must be implemented by the child)
Definition: Monomial.cpp:22
virtual double Normalization(unsigned int polyOrder) const
virtual double BasisEvaluate(int const order, double const x) const override
Evaluate the specific polynomial type (must be implemented by the child)
static std::shared_ptr< OrthogonalPolynomial > Construct(std::string const &polyName)
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
virtual double Normalization(unsigned int polyOrder) const override
virtual double DerivativeEvaluate(int const polyOrder, int const derivOrder, double const x) const override
virtual double Normalization(unsigned int polyOrder) const override
void PolynomialsWrapper(pybind11::module &m)