version 0.4.1
python/dirichletvalues/dirichletvalues.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 <cstdlib>
12#include <string>
13
14#include "dune/common/classname.hh"
15#include <dune/functions/functionspacebases/lagrangebasis.hh>
16#include <dune/functions/functionspacebases/powerbasis.hh>
17#include <dune/grid/yaspgrid.hh>
18#include <dune/python/common/typeregistry.hh>
19#include <dune/python/functions/globalbasis.hh>
20#include <dune/python/functions/subspacebasis.hh>
21#include <dune/python/pybind11/eigen.h>
22#include <dune/python/pybind11/functional.h>
23#include <dune/python/pybind11/pybind11.h>
24#include <dune/python/pybind11/stl.h>
25#include <dune/python/pybind11/stl_bind.h>
26
28
29namespace Ikarus::Python {
30
31namespace Impl {
32 using FixBoundaryDOFsWithGlobalIndexFunction = std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int)>;
33
34 template <typename LV>
35 using FixBoundaryDOFsWithLocalViewFunction = std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LV&)>;
36
37 template <typename LV, typename IS>
38 using FixBoundaryDOFsWithIntersectionFunction =
39 std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LV&, const IS&)>;
40
41 template <typename Basis>
42 auto registerSubSpaceLocalView() {
43 pybind11::module scopedf = pybind11::module::import("dune.functions");
44 using LocalViewWrapper = Dune::Python::LocalViewWrapper<Basis>;
45
46 auto includes = Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"};
47
48 Dune::Python::insertClass<Basis>(scopedf, "SubspaceBasis_" + Dune::className<typename Basis::PrefixPath>(),
49 Dune::Python::GenerateTypeName(Dune::className<Basis>()), includes);
50
51 auto [lv, isNew] = Dune::Python::insertClass<LocalViewWrapper>(
52 scopedf, "LocalView_" + Dune::className<typename Basis::PrefixPath>(),
53 Dune::Python::GenerateTypeName("Dune::Python::LocalViewWrapper", Dune::MetaType<Basis>()), includes);
54 if (isNew) {
55 lv.def("bind", &LocalViewWrapper::bind);
56 lv.def("unbind", &LocalViewWrapper::unbind);
57 lv.def("index", [](const LocalViewWrapper& localView, int index) { return localView.index(index); });
58 lv.def("__len__", [](LocalViewWrapper& self) -> int { return self.size(); });
59
60 Dune::Python::Functions::registerTree<typename LocalViewWrapper::Tree>(lv);
61 lv.def("tree", [](const LocalViewWrapper& view) { return view.tree(); });
62 }
63 }
64} // namespace Impl
65
66template <class DirichletValues>
67void forwardCorrectFunction(DirichletValues& dirichletValues, const pybind11::function& functor, auto&& cppFunction) {
68 using Basis = typename DirichletValues::Basis;
69 using Intersection = typename Basis::GridView::Intersection;
70 using BackendType = typename DirichletValues::BackendType;
71 using MultiIndex = typename Basis::MultiIndex;
72
73 // Disambiguate by number of arguments
74 pybind11::module inspect_module = pybind11::module::import("inspect");
75 pybind11::object result = inspect_module.attr("signature")(functor).attr("parameters");
76 size_t numParams = pybind11::len(result);
77
78 if (numParams == 2) {
79 auto function = functor.template cast<const Impl::FixBoundaryDOFsWithGlobalIndexFunction>();
80 auto lambda = [&](BackendType& vec, const MultiIndex& indexGlobal) { function(vec.vector(), indexGlobal); };
81 cppFunction(lambda);
82
83 } else if (numParams == 3) {
84 auto lambda = [&](BackendType& vec, int localIndex, auto&& lv) {
85 using SubSpaceBasis = typename std::remove_cvref_t<decltype(lv)>::GlobalBasis;
86 Impl::registerSubSpaceLocalView<SubSpaceBasis>();
87
88 using SubSpaceLocalViewWrapper = Dune::Python::LocalViewWrapper<SubSpaceBasis>;
89 auto lvWrapper = SubSpaceLocalViewWrapper(lv);
90
91 auto function =
92 functor.template cast<const Impl::FixBoundaryDOFsWithLocalViewFunction<SubSpaceLocalViewWrapper>>();
93 function(vec.vector(), localIndex, lvWrapper);
94 };
95 cppFunction(lambda);
96
97 } else if (numParams == 4) {
98 auto lambda = [&](BackendType& vec, int localIndex, auto&& lv, const Intersection& intersection) {
99 using SubSpaceBasis = typename std::remove_cvref_t<decltype(lv)>::GlobalBasis;
100 Impl::registerSubSpaceLocalView<SubSpaceBasis>();
101
102 using SubSpaceLocalViewWrapper = Dune::Python::LocalViewWrapper<SubSpaceBasis>;
103 auto lvWrapper = SubSpaceLocalViewWrapper(lv);
104
105 auto function = functor.template cast<
106 const Impl::FixBoundaryDOFsWithIntersectionFunction<SubSpaceLocalViewWrapper, Intersection>>();
107 function(vec.vector(), localIndex, lvWrapper, intersection);
108 };
109 cppFunction(lambda);
110
111 } else {
112 DUNE_THROW(Dune::NotImplemented, "fixBoundaryDOFs: A function with this signature is not supported");
113 }
114}
115
146template <class DirichletValues, class... options>
147void registerDirichletValues(pybind11::handle scope, pybind11::class_<DirichletValues, options...> cls) {
148 using pybind11::operator""_a;
149
150 using Basis = typename DirichletValues::Basis;
151 using BackendType = typename DirichletValues::BackendType;
152 using FlagsType = typename DirichletValues::FlagsType;
153 using MultiIndex = typename Basis::MultiIndex;
154 using LocalView = typename Basis::LocalView;
155 using Intersection = typename Basis::GridView::Intersection;
156
157 pybind11::module scopedf = pybind11::module::import("dune.functions");
158 using LocalViewWrapper = Dune::Python::LocalViewWrapper<Basis>;
159
160 auto includes = Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"};
161 auto [lv, isNew] = Dune::Python::insertClass<LocalViewWrapper>(
162 scopedf, "LocalView", Dune::Python::GenerateTypeName("Dune::Python::LocalViewWrapper", Dune::MetaType<Basis>()),
163 includes);
164
165 if (isNew) {
166 lv.def("bind", &LocalViewWrapper::bind);
167 lv.def("unbind", &LocalViewWrapper::unbind);
168 lv.def("index", [](const LocalViewWrapper& localView, int index) { return localView.index(index); });
169 lv.def("__len__", [](LocalViewWrapper& self) -> int { return self.size(); });
170
171 Dune::Python::Functions::registerTree<typename LocalViewWrapper::Tree>(lv);
172 lv.def("tree", [](const LocalViewWrapper& view) { return view.tree(); });
173 }
174
175 cls.def(pybind11::init([](const Basis& basis) { return new DirichletValues(basis); }), pybind11::keep_alive<1, 2>());
176
177 cls.def_property_readonly("container", &DirichletValues::container);
178 cls.def_property_readonly("size", &DirichletValues::size);
179 cls.def("__len__", &DirichletValues::size);
180 cls.def_property_readonly("fixedDOFsize", &DirichletValues::fixedDOFsize);
181 cls.def("isConstrained", [](DirichletValues& self, std::size_t i) -> bool { return self.isConstrained(i); });
182 cls.def("setSingleDOF", [](DirichletValues& self, std::size_t i, bool flag) { self.setSingleDOF(i, flag); });
183 cls.def("isConstrained", [](DirichletValues& self, MultiIndex i) -> bool { return self.isConstrained(i); });
184 cls.def("setSingleDOF", [](DirichletValues& self, MultiIndex i, bool flag) { self.setSingleDOF(i, flag); });
185 cls.def("reset", &DirichletValues::reset);
186
187 cls.def("fixDOFs",
188 [](DirichletValues& self, const std::function<void(const Basis&, Eigen::Ref<Eigen::VectorX<bool>>)>& f) {
189 auto lambda = [&](const Basis& basis, BackendType& vec) { f(basis, vec.vector()); };
190 self.fixDOFs(lambda);
191 });
192}
193
194} // namespace Ikarus::Python
Definition of the LinearElastic class for finite element mechanics computations.
void registerDirichletValues(pybind11::handle scope, pybind11::class_< DirichletValues, options... > cls)
Register Python bindings for a DirichletValues class.
Definition: python/dirichletvalues/dirichletvalues.hh:147
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:82
Definition: flatassembler.hh:21
void forwardCorrectFunction(DirichletValues &dirichletValues, const pybind11::function &functor, auto &&cppFunction)
Definition: python/dirichletvalues/dirichletvalues.hh:67
def basis(gv, tree)
Definition: basis.py:10
def dirichletValues(basis)
Definition: dirichlet_values.py:38
Class for handling Dirichlet boundary conditions in Ikarus.
Definition: utils/dirichletvalues.hh:67
auto size() const
Definition: utils/dirichletvalues.hh:181
void setSingleDOF(const MultiIndex i, bool flag)
Fixes and unfixes (set boolean value to true or false) a specific degree of freedom.
Definition: utils/dirichletvalues.hh:135
bool isConstrained(const MultiIndex &multiIndex) const
Definition: utils/dirichletvalues.hh:166
std::remove_cvref_t< B > Basis
Definition: utils/dirichletvalues.hh:69
decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >())) BackendType
Definition: utils/dirichletvalues.hh:72
void fixDOFs(F &&f)
Function to fix (set boolean values to true or false) degrees of freedom.
Definition: utils/dirichletvalues.hh:123
FC FlagsType
Definition: utils/dirichletvalues.hh:70
void reset()
Resets all degrees of freedom.
Definition: utils/dirichletvalues.hh:155
auto fixedDOFsize() const
Definition: utils/dirichletvalues.hh:178
auto & container() const
Definition: utils/dirichletvalues.hh:184