3#ifndef DUNE_PYTHON_LOCALFUNCTIONS_LOCALFINITEELEMENT_HH
4#define DUNE_PYTHON_LOCALFUNCTIONS_LOCALFINITEELEMENT_HH
6#include <dune/python/pybind11/pybind11.h>
8#include <dune/common/visibility.hh>
17template<
typename LocalBasis>
18DUNE_EXPORT
auto registerLocalBasis(pybind11::handle scope)
20 static auto cls = pybind11::class_<LocalBasis>(scope,
"LocalBasis");
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);
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);
39DUNE_EXPORT
auto registerLocalKey(pybind11::handle scope)
41 static auto cls = pybind11::class_<LocalKey>(scope,
"LocalKey");
45 cls.def_property(
"index",
48 cls.def(
"__lt__", &LocalKey::operator<);
55template<
typename LocalFiniteElement>
58 static auto cls = pybind11::class_<LocalFiniteElement>(scope, name);
60 detail::registerLocalBasis<typename LocalFiniteElement::Traits::LocalBasisType>(cls);
Definition bdfmcube.hh:18
Definition python/localfunctions/localfiniteelement.hh:13
DUNE_EXPORT auto registerLocalFiniteElement(pybind11::handle scope, const char *name="LocalFiniteElement")
Definition python/localfunctions/localfiniteelement.hh:56
Type erasure class storing a local finite element.
Definition localfunctions/common/localfiniteelement.hh:40
const Traits::LocalBasisType & localBasis() const
Access the LocalBasis of the stored local finite element.
Definition localfunctions/common/localfiniteelement.hh:156
unsigned int size() const
Get the number of basis functions of the stored local finite element.
Definition localfunctions/common/localfiniteelement.hh:188
const GeometryType & type() const
Get the GeometryType of the stored local finite element.
Definition localfunctions/common/localfiniteelement.hh:198
Describe position of one degree of freedom.
Definition localkey.hh:24
constexpr unsigned int index() const noexcept
Return offset within subentity.
Definition localkey.hh:70
constexpr unsigned int codim() const noexcept
Return codim of associated entity.
Definition localkey.hh:63
constexpr unsigned int subEntity() const noexcept
Return number of associated subentity.
Definition localkey.hh:56