Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F91370581
lm_python_domains.hh
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sun, Nov 10, 10:47
Size
3 KB
Mime Type
text/x-c++
Expires
Tue, Nov 12, 10:47 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
22252030
Attached To
rLIBMULTISCALE LibMultiScale
lm_python_domains.hh
View Options
#ifndef __LM_PYTHON_DOMAINS_HH__
#define __LM_PYTHON_DOMAINS_HH__
/* -------------------------------------------------------------------------- */
#include <pybind11/eigen.h>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
/* -------------------------------------------------------------------------- */
#include "domain_multiscale.hh"
/* -------------------------------------------------------------------------- */
namespace py = pybind11;
__BEGIN_LIBMULTISCALE__
/* -------------------------------------------------------------------------- */
template <typename Class>
inline void declare_domain(py::module &m, const char *class_name) {
py::class_<Class, DomainInterface, std::shared_ptr<Class>>(
m, class_name, py::multiple_inheritance())
.def(py::init<LMID, CommGroup &>())
.def("init",
[](Class &self) {
self.checkAllKeywordsAreParsed();
self.init();
})
.def("getContainer", &Class::getContainer,
py::return_value_policy::reference)
.def("updateGradient", &Class::updateGradient);
std::string getter_name = std::string(class_name) + "Container";
using container_class =
std::decay_t<decltype(std::declval<Class>().getContainer())>;
py::class_<container_class, ContainerInterface,
std::shared_ptr<container_class>>(m, getter_name.c_str());
}
/* -------------------------------------------------------------------------- */
inline void declare_domains(py::module &m) {
auto clean_class_name = [](const std::string &arg) {
std::string res = arg;
res.erase(std::remove(res.begin(), res.end(), '<'), res.end());
res.erase(std::remove(res.begin(), res.end(), '>'), res.end());
return res;
};
py::class_<DomainMultiScale>(m, "DomainMultiScale")
.def_static("getManager", &DomainMultiScale::getManager,
py::return_value_policy::reference)
.def("build", &DomainMultiScale::build)
.def("getObject", &DomainMultiScale::getObject,
py::return_value_policy::reference)
.def("coupling", &DomainMultiScale::coupling)
.def("destroy", &DomainMultiScale::destroy)
.def("addObject", py::overload_cast<std::shared_ptr<DomainInterface>>(
&DomainMultiScale::addObject));
py::class_<DomainInterface, Parsable, std::shared_ptr<DomainInterface>>(
m, "DomainInterface", py::multiple_inheritance());
#define PYTHON_DOMAIN(n, data, obj) \
{ \
using _class = BOOST_PP_TUPLE_ELEM(3, 0, obj); \
std::string domain_name = \
clean_class_name(BOOST_PP_STRINGIZE(BOOST_PP_TUPLE_ELEM(3, 0, obj))); \
declare_domain<_class>(m, domain_name.c_str()); \
}
BOOST_PP_SEQ_FOR_EACH(PYTHON_DOMAIN, f, LIST_ATOM_MODEL);
BOOST_PP_SEQ_FOR_EACH(PYTHON_DOMAIN, f, LIST_DD_MODEL);
BOOST_PP_SEQ_FOR_EACH(PYTHON_DOMAIN, f, LIST_CONTINUUM_MODEL);
#undef PYTHON_DOMAIN
}
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
#endif //__LM_PYTHON_DOMAINS_HH__
Event Timeline
Log In to Comment