dune-localfunctions  2.8.0
python/localfunctions/localfiniteelement.hh
Go to the documentation of this file.
1 #ifndef DUNE_PYTHON_LOCALFUNCTIONS_LOCALFINITEELEMENT_HH
2 #define DUNE_PYTHON_LOCALFUNCTIONS_LOCALFINITEELEMENT_HH
3 
4 #include <dune/python/pybind11/pybind11.h>
5 
6 #include <dune/common/visibility.hh>
9 
10 namespace Dune {
11 namespace Python {
12 
13 namespace detail {
14 
15 template<typename LocalBasis>
16 DUNE_EXPORT auto registerLocalBasis(pybind11::handle scope)
17 {
18  static auto cls = pybind11::class_<LocalBasis>(scope, "LocalBasis");
19 
20  cls.def("__len__", [](const LocalBasis& basis) { return basis.size(); });
21  cls.def_property_readonly("order", [](const LocalBasis& basis) { return basis.order(); });
22  cls.def("evaluateFunction",
23  [](const LocalBasis& basis, const typename LocalBasis::Traits::DomainType& in) {
24  std::vector<typename LocalBasis::Traits::RangeType> out;
25  basis.evaluateFunction(in, out);
26  return out;
27  });
28  cls.def("evaluateJacobian",
29  [](const LocalBasis& basis, const typename LocalBasis::Traits::DomainType& in) {
30  std::vector<typename LocalBasis::Traits::JacobianType> out;
31  basis.evaluateJacobian(in, out);
32  return out;
33  });
34  return cls;
35 }
36 
37 DUNE_EXPORT auto registerLocalKey(pybind11::handle scope)
38 {
39  static auto cls = pybind11::class_<LocalKey>(scope, "LocalKey");
40 
41  cls.def_property_readonly("subEntity", &LocalKey::subEntity);
42  cls.def_property_readonly("codim", &LocalKey::codim);
43  cls.def_property("index",
44  [](const LocalKey& key) { return key.index(); },
45  [](LocalKey& key, unsigned int index) { key.index(index); });
46  cls.def("__lt__", &LocalKey::operator<);
47 
48  return cls;
49 }
50 
51 } /* namespace detail */
52 
53 template<typename LocalFiniteElement>
54 DUNE_EXPORT auto registerLocalFiniteElement(pybind11::handle scope, const char* name = "LocalFiniteElement")
55 {
56  static auto cls = pybind11::class_<LocalFiniteElement>(scope, name);
57 
58  detail::registerLocalBasis<typename LocalFiniteElement::Traits::LocalBasisType>(cls);
59 
60  cls.def_property_readonly("localBasis", &LocalFiniteElement::localBasis, pybind11::return_value_policy::reference_internal);
61  // cls.def_property_readonly("localCoefficients", &LocalFiniteElement::localCoefficients);
62  // cls.def_property_readonly("localInterpolation", &LocalFiniteElement::localInterpolation);
63  cls.def("__len__", &LocalFiniteElement::size);
64  cls.def_property_readonly("type", &LocalFiniteElement::type);
65 
66  return cls;
67 }
68 
69 
70 } /* namespace Python */
71 } /* namespace Dune */
72 
73 #endif
Definition: bdfmcube.hh:16
DUNE_EXPORT auto registerLocalFiniteElement(pybind11::handle scope, const char *name="LocalFiniteElement")
Definition: python/localfunctions/localfiniteelement.hh:54
Describe position of one degree of freedom.
Definition: localkey.hh:21
unsigned int index() const
Return offset within subentity.
Definition: localkey.hh:66
unsigned int codim() const
Return codim of associated entity.
Definition: localkey.hh:60
unsigned int subEntity() const
Return number of associated subentity.
Definition: localkey.hh:54