MUQ  0.4.3
MultiIndicesWrapper.cpp
Go to the documentation of this file.
1 #include "AllClassWrappers.h"
2 
7 
8 #include <pybind11/pybind11.h>
9 #include <pybind11/stl.h>
10 #include <pybind11/eigen.h>
11 
12 #include <string>
13 
14 #include <functional>
15 #include <vector>
16 
17 using namespace muq::Utilities::PythonBindings;
18 namespace py = pybind11;
19 
21 {
22  py::class_<MultiIndex, std::shared_ptr<MultiIndex>> multiI(m, "MultiIndex");
23  multiI
24  .def(py::init<unsigned>())
25  .def(py::init<Eigen::RowVectorXi const&>())
26  .def(py::init<std::initializer_list<unsigned> const&>())
27  .def_static("Copy", &MultiIndex::Copy)
28  .def("GetVector", &MultiIndex::GetVector)
29  .def("Sum", &MultiIndex::Sum)
30  .def("Max", &MultiIndex::Max)
31  .def("SetValue", &MultiIndex::SetValue)
32  .def("GetValue", &MultiIndex::GetValue)
33  .def("SetLength", &MultiIndex::SetLength)
34  .def("GetLength", &MultiIndex::GetLength)
35  .def("GetNzBegin", &MultiIndex::GetNzBegin)
36  .def("GetNzEnd", &MultiIndex::GetNzEnd);
37 
38  py::class_<MultiIndexLimiter, std::shared_ptr<MultiIndexLimiter>> multiILim(m, "MultiIndexLimiter");
39  multiILim
40  .def("IsFeasible", &MultiIndexLimiter::IsFeasible);
41 
42  py::class_<TotalOrderLimiter, MultiIndexLimiter, std::shared_ptr<TotalOrderLimiter>> totalOrd(m, "TotalOrderLimiter");
43  totalOrd
44  .def(py::init<unsigned int>());
45 
46  py::class_<DimensionLimiter, MultiIndexLimiter, std::shared_ptr<DimensionLimiter>> dimLim(m, "DimensionLimiter");
47  dimLim
48  .def(py::init<unsigned int, unsigned int>())
49  .def("IsFeasible", &DimensionLimiter::IsFeasible);
50 
51  py::class_<MaxOrderLimiter, MultiIndexLimiter, std::shared_ptr<MaxOrderLimiter>> maxLim(m, "MaxOrderLimiter");
52  maxLim
53  .def(py::init<unsigned int>())
54  .def(py::init<Eigen::VectorXi const&>())
55  .def("IsFeasible", &MaxOrderLimiter::IsFeasible);
56 
57  py::class_<NoLimiter, MultiIndexLimiter, std::shared_ptr<NoLimiter>> noLim(m, "NoLimiter");
58  noLim
59  .def(py::init<>())
60  .def("IsFeasible", &NoLimiter::IsFeasible);
61 
62  py::class_<AndLimiter, MultiIndexLimiter, std::shared_ptr<AndLimiter>> andLim(m, "AndLimiter");
63  andLim
64  .def(py::init<std::shared_ptr<MultiIndexLimiter>, std::shared_ptr<MultiIndexLimiter>>())
65  .def("IsFeasible", &AndLimiter::IsFeasible);
66 
67  py::class_<OrLimiter, MultiIndexLimiter, std::shared_ptr<OrLimiter>> orLim(m, "OrLimiter");
68  orLim
69  .def(py::init<std::shared_ptr<MultiIndexLimiter>, std::shared_ptr<MultiIndexLimiter>>())
70  .def("IsFeasible", &OrLimiter::IsFeasible);
71 
72  py::class_<XorLimiter, MultiIndexLimiter, std::shared_ptr<XorLimiter>> xorLim(m, "XorLimiter");
73  xorLim
74  .def(py::init<std::shared_ptr<MultiIndexLimiter>, std::shared_ptr<MultiIndexLimiter>>())
75  .def("IsFeasible", &XorLimiter::IsFeasible);
76 
77  py::class_<MultiIndexFactory, std::shared_ptr<MultiIndexFactory>> multiIFac(m, "MultiIndexFactory");
78  multiIFac
79  .def_static("CreateTotalOrder", &MultiIndexFactory::CreateTotalOrder, py::arg("length"), py::arg("maxOrder"), py::arg("minOrder")=0, py::arg("limiter")=std::make_shared<NoLimiter>())
80  .def_static("CreateTriTotalOrder", &MultiIndexFactory::CreateTriTotalOrder, py::arg("length"), py::arg("maxOrder"), py::arg("minOrder")=0, py::arg("limiter")=std::make_shared<NoLimiter>())
81  .def_static("CreateHyperbolic", &MultiIndexFactory::CreateHyperbolic, py::arg("length"), py::arg("maxOrder"), py::arg("q")=0, py::arg("limiter")=std::make_shared<NoLimiter>())
82  .def_static("CreateTriHyperbolic", &MultiIndexFactory::CreateTriHyperbolic, py::arg("length"), py::arg("maxOrder"), py::arg("q")=0, py::arg("limiter")=std::make_shared<NoLimiter>())
83  .def_static("CreateFullTensor", (std::shared_ptr<MultiIndexSet>(*)(unsigned int const, unsigned int const,std::shared_ptr<MultiIndexLimiter>)) &MultiIndexFactory::CreateFullTensor, py::arg("length"), py::arg("order"), py::arg("limiter")=std::make_shared<NoLimiter>())
84  .def_static("CreateFullTensor", (std::shared_ptr<MultiIndexSet>(*)(const Eigen::RowVectorXi&, std::shared_ptr<MultiIndexLimiter>)) &MultiIndexFactory::CreateFullTensor, py::arg("orders"), py::arg("limiter")=std::make_shared<NoLimiter>())
85  .def_static("CreateAnisotropic", (std::shared_ptr<MultiIndexSet>(*)(const Eigen::RowVectorXf&, const double)) &MultiIndexFactory::CreateAnisotropic, py::arg("weights"), py::arg("epsilon"))
86  //.def("CentralMoment", &SampleCollection::CentralMoment, py::arg("order"), py::arg("blockDim") = -1)
87  .def_static("CreateTriHyperbolic", &MultiIndexFactory::CreateTriHyperbolic);
88 
89  py::class_<MultiIndexSet,std::shared_ptr<MultiIndexSet>> multiSet(m, "MultiIndexSet");
90  multiSet
91  .def(py::init<const unsigned , std::shared_ptr<MultiIndexLimiter>>())
92  .def_static("CloneExisting", &MultiIndexSet::CloneExisting)
93  .def("SetLimiter", &MultiIndexSet::SetLimiter)
94  .def("GetLimiter", &MultiIndexSet::GetLimiter)
95  .def("IndexToMulti", &MultiIndexSet::IndexToMulti)
96  .def("MultiToIndex", &MultiIndexSet::MultiToIndex)
97  .def("GetMultiLength", &MultiIndexSet::GetMultiLength)
98  .def("GetAllMultiIndices", &MultiIndexSet::GetAllMultiIndices)
99  .def("GetMaxOrders", &MultiIndexSet::GetMaxOrders)
100  .def("at", &MultiIndexSet::at)
101  .def("Size", &MultiIndexSet::Size)
102  .def("Union", &MultiIndexSet::Union)
103  .def("Activate", (void (MultiIndexSet::*)(std::shared_ptr<MultiIndex> const&)) &MultiIndexSet::Activate)
104  .def("AddActive", &MultiIndexSet::AddActive)
105  .def("Expand", &MultiIndexSet::Expand)
106  .def("ForciblyExpand", &MultiIndexSet::ForciblyExpand)
107  .def("ForciblyActivate", (std::vector<unsigned> (MultiIndexSet::*)(std::shared_ptr<MultiIndex> const&)) &MultiIndexSet::ForciblyActivate)
108  .def("GetAdmissibleForwardNeighbors", &MultiIndexSet::GetAdmissibleForwardNeighbors)
109  .def("IsAdmissible", (bool (MultiIndexSet::*)(std::shared_ptr<MultiIndex> const&) const) &MultiIndexSet::IsAdmissible)
110  .def("IsExpandable", &MultiIndexSet::IsExpandable)
111  .def("IsActive", (bool (MultiIndexSet::*)(std::shared_ptr<MultiIndex> const&) const) &MultiIndexSet::IsActive)
112  .def("ToHDF5", (void (MultiIndexSet::*)(std::string, std::string) const) &MultiIndexSet::ToHDF5, py::arg("filename"),py::arg("dsetName")="/multiindices")
113  .def_static("FromHDF5", [](std::string filename, std::string dsetName){return MultiIndexSet::FromHDF5(filename, dsetName);})
114  .def_static("FromHDF5", [](std::string filename){return MultiIndexSet::FromHDF5(filename, "/multiindices");});
115 
116 // .def_static("FromHDF5", (std::shared_ptr<MultiIndexSet> (MultiIndexSet::*)(std::string, std::string)) &MultiIndexSet::FromHDF5, py::arg("filename"), py::arg("dsetName")=std::string("/multiindices"));
117 
118 
119 }
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
static std::shared_ptr< MultiIndexSet > CreateTotalOrder(unsigned int const length, unsigned int const maxOrder, unsigned int const minOrder=0, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a total order limited MultiIndex.
static std::vector< std::shared_ptr< MultiIndexSet > > CreateTriHyperbolic(unsigned int const length, unsigned int const maxOrder, const double q=0.5, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
static std::shared_ptr< MultiIndexSet > CreateFullTensor(unsigned int const length, unsigned int const order, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a full tensor product multiindex set.
static std::shared_ptr< MultiIndexSet > CreateHyperbolic(unsigned int const length, unsigned int const maxOrder, const double q=0.5, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a general hyperbolic MultiIndex.
static std::vector< std::shared_ptr< MultiIndexSet > > CreateTriTotalOrder(unsigned int const length, unsigned int const maxOrder, unsigned int const minOrder=0, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
static std::shared_ptr< MultiIndexSet > CreateAnisotropic(const Eigen::RowVectorXf &weights, const double epsilon)
Construct an anisotropic multiindex set based on a priori information on the importance of each dimen...
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const =0
virtual bool IsActive(std::shared_ptr< MultiIndex > const &multiIndex) const
Return true if the multiIndex is active.
virtual bool IsExpandable(unsigned int activeIndex) const
Return true if one of the forward neighbors of index is admissible but not active.
virtual std::vector< unsigned > Expand(unsigned int activeIndex)
virtual unsigned int Size() const
static std::shared_ptr< MultiIndexSet > FromHDF5(std::string filename, std::string dsetName="/multiiindices")
Loads a multiindex set from an HDF5 file.
static std::shared_ptr< MultiIndexSet > CloneExisting(std::shared_ptr< MultiIndexSet > const &original)
NOTE: does not perform a deep copy of the multiindices themselves, only pointers to the multiindices.
virtual std::vector< std::shared_ptr< MultiIndex > > GetAllMultiIndices() const
Definition: MultiIndexSet.h:97
virtual int Union(const MultiIndexSet &rhs)
Add all terms in rhs to this instance.
virtual void SetLimiter(std::shared_ptr< MultiIndexLimiter > const &limiterIn)
virtual Eigen::VectorXi GetMaxOrders() const
virtual unsigned int GetMultiLength() const
virtual std::vector< unsigned > ForciblyActivate(std::shared_ptr< MultiIndex > const &multiIndex)
virtual int MultiToIndex(std::shared_ptr< MultiIndex > const &input) const
virtual void Activate(std::shared_ptr< MultiIndex > const &multiIndex)
void ToHDF5(std::string filename, std::string dsetName="/multiindices") const
Saves the multiindex set to a group in an HDF5 file.
virtual bool IsAdmissible(std::shared_ptr< MultiIndex > const &multiIndex) const
Determines whether the input multiIndex is currently admissible.
virtual std::shared_ptr< MultiIndex > const & at(int activeIndex)
virtual std::vector< unsigned > ForciblyExpand(unsigned int const activeIndex)
virtual int AddActive(std::shared_ptr< MultiIndex > const &newNode)
virtual std::shared_ptr< MultiIndex > const & IndexToMulti(unsigned activeIndex) const
virtual std::vector< std::shared_ptr< MultiIndex > > GetAdmissibleForwardNeighbors(unsigned int activeIndex)
virtual std::shared_ptr< MultiIndexLimiter > GetLimiter() const
Definition: MultiIndexSet.h:88
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
virtual bool IsFeasible(std::shared_ptr< MultiIndex > multi) const override
void MultiIndicesWrapper(pybind11::module &m)