version 0.4.1
fe.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10
11#include <dune/python/common/typeregistry.hh>
12#include <dune/python/functions/globalbasis.hh>
13#include <dune/python/pybind11/eigen.h>
14#include <dune/python/pybind11/functional.h>
15#include <dune/python/pybind11/pybind11.h>
16#include <dune/python/pybind11/stl.h>
17
20
21namespace Ikarus::Python {
22
50template <class FE, class... options>
51void registerCalculateAt(pybind11::handle scope, pybind11::class_<FE, options...> cls, auto resultTypesTuple) {
52 using Traits = typename FE::Traits;
53 using FERequirements = typename FE::Requirement;
54 cls.def(
55 "calculateAt",
56 [&](FE& self, const FERequirements& req, const Dune::FieldVector<double, Traits::mydim>& local,
57 std::string resType) {
58 Eigen::VectorXd result;
59 bool success = false;
60 Dune::Hybrid::forEach(resultTypesTuple, [&]<typename RT>(RT i) {
61 if (resType == toString(i)) {
62 success = true;
63 result = self.template calculateAt<RT::template Rebind>(req, local).asVec();
64 }
65 });
66 if (success)
67 return result;
68 DUNE_THROW(Dune::NotImplemented, "Element " + Dune::className<FE>() + " doesn't support ResultType " + resType);
69 },
70 pybind11::arg("feRequirements"), pybind11::arg("local"), pybind11::arg("resultType"));
71}
72
86template <class FE, class... options>
87void registerFE(pybind11::handle scope, pybind11::class_<FE, options...> cls) {
88 using BH = typename FE::BasisHandler;
89 using GridElement = typename FE::GridElement;
90 using Requirement = typename FE::Requirement;
91 using FlatBasis = typename FE::Traits::FlatBasis;
92
93 int index = 0;
94 cls.def(pybind11::init([](const BH& basisHandler, typename FE::PreTuple argsTuple) {
95 auto unpackTuple = [&]<typename... Arg>(Arg&&... args) {
96 return new FE(basisHandler, std::forward<Arg>(args)...);
97 };
98 return std::apply(unpackTuple, argsTuple);
99 }),
100 pybind11::keep_alive<1, 2>());
101
102 cls.def("bind", [](FE& self, const GridElement& e) { self.bind(e); });
103 cls.def("calculateScalar", [](FE& self, const Requirement& req, Ikarus::ScalarAffordance affordance) {
104 return calculateScalar(self, req, affordance);
105 });
106 cls.def("calculateVector", [](FE& self, const Requirement& req, Ikarus::VectorAffordance affordance,
107 Eigen::Ref<Eigen::VectorXd> vec) { calculateVector(self, req, affordance, vec); });
108 cls.def(
109 "calculateMatrix",
110 [](FE& self, const Requirement& req, Ikarus::MatrixAffordance affordance, Eigen::Ref<Eigen::MatrixXd> mat) {
111 calculateMatrix(self, req, affordance, mat);
112 },
113 pybind11::arg("Requirement"), pybind11::arg("MatrixAffordance"), pybind11::arg("elementMatrix").noconvert());
114
115 pybind11::module scopedf = pybind11::module::import("dune.functions");
116 using LocalViewWrapper = Dune::Python::LocalViewWrapper<FlatBasis>;
117
118 auto includes = Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"};
119 auto [lv, isNew] = Dune::Python::insertClass<LocalViewWrapper>(
120 scopedf, "LocalViewWrapper",
121 Dune::Python::GenerateTypeName("Dune::Python::LocalViewWrapperWrapper", Dune::MetaType<FlatBasis>()), includes);
122
123 if (isNew) {
124 lv.def("bind", &LocalViewWrapper::bind);
125 lv.def("unbind", &LocalViewWrapper::unbind);
126 lv.def("index", [](const LocalViewWrapper& localView, int index) { return localView.index(index); });
127 lv.def("__len__", [](LocalViewWrapper& self) -> int { return self.size(); });
128
129 Dune::Python::Functions::registerTree<typename LocalViewWrapper::Tree>(lv);
130 lv.def("tree", [](const LocalViewWrapper& view) { return view.tree(); });
131 }
132
133 cls.def(
134 "localView",
135 [](FE& self) -> LocalViewWrapper {
136 auto lvWrapped = LocalViewWrapper(self.localView().globalBasis());
137 // this can be simplified when https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418
138 // becomes available
139 // pybind11::object obj = pybind11::cast(self.localView().element());
140 lvWrapped.base().bind(self.localView().element());
141 return lvWrapped;
142 },
143 pybind11::keep_alive<0, 1>());
144
145 registerCalculateAt(scope, cls, typename FE::SupportedResultTypes());
146 registerFERequirement(scope, cls);
147}
148
149} // namespace Ikarus::Python
Contains the definition of the FEFactory class.
void registerFE(pybind11::handle scope, pybind11::class_< FE, options... > cls)
Register Python bindings for the FE class.
Definition: fe.hh:87
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:82
constexpr std::string toString(DBCOption _e)
Definition: dirichletbcenforcement.hh:7
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
Definition: flatassembler.hh:21
void registerCalculateAt(pybind11::handle scope, pybind11::class_< FE, options... > cls, auto resultTypesTuple)
Registers the calculateAt method for a finite element class in Python.
Definition: fe.hh:51
void registerFERequirement(pybind11::handle scope, pybind11::class_< FE, options... > cls)
Definition: registerferequirements.hh:14
FE class is a base class for all finite elements.
Definition: febase.hh:79
const LocalView & localView() const
Get the const reference to the local view.
Definition: febase.hh:134
Traits::BasisHandler BasisHandler
Type of the basisHandler.
Definition: febase.hh:88
std::tuple< typename Skills< PreFE, FE >::Pre... > PreTuple
Definition: febase.hh:97
void bind(const GridElement &element)
Convenient function to bind the local view to the element.
Definition: febase.hh:113
PreFE::Traits Traits
Type traits.
Definition: febase.hh:87
typename Traits::Element GridElement
Type of the grid element.
Definition: febase.hh:92
Class representing the requirements for finite element calculations.
Definition: ferequirements.hh:252
typename BasisHandler::FlatBasis FlatBasis
Type of the flat basis.
Definition: fetraits.hh:33
RequirementType< requirementDetected >::type Requirement
Definition: mixin.hh:89
decltype(std::tuple_cat(computeSupportedResultTypes< Skills...< PreFE, typename PreFE::template FE< Skills... ... > > >()...)) SupportedResultTypes
Type alias for the supported result types by the mixin.
Definition: mixin.hh:64
Definition: utils/dirichletvalues.hh:32