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