Dune Core Modules (2.9.0)

localfiniteelement.hh
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
9#include <dune/localfunctions/common/localkey.hh>
10#include <dune/localfunctions/common/virtualinterface.hh>
11
12namespace Dune {
13namespace Python {
14
15namespace detail {
16
17template<typename LocalBasis>
18DUNE_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
39DUNE_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
55template<typename LocalFiniteElement>
56DUNE_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
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
Dune namespace.
Definition: alignedallocator.hh:13
Definition of macros controlling symbol visibility at the ABI level.
#define DUNE_EXPORT
Export a symbol as part of the public ABI.
Definition: visibility.hh:20
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jan 6, 23:30, 2025)