version 0.4
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 <dune/functions/functionspacebases/lagrangebasis.hh>
12#include <dune/functions/functionspacebases/powerbasis.hh>
13#include <dune/grid/yaspgrid.hh>
14#include <dune/python/common/typeregistry.hh>
15#include <dune/python/functions/globalbasis.hh>
16#include <dune/python/pybind11/eigen.h>
17#include <dune/python/pybind11/functional.h>
18#include <dune/python/pybind11/pybind11.h>
19#include <dune/python/pybind11/stl.h>
20#include <dune/python/pybind11/stl_bind.h>
21
23// PYBIND11_MAKE_OPAQUE(std::vector<bool>);
24namespace Ikarus::Python {
25
48 template <class DirichletValues, class... options>
49 void registerDirichletValues(pybind11::handle scope, pybind11::class_<DirichletValues, options...> cls) {
50 using pybind11::operator""_a;
51
52 using Basis = typename DirichletValues::Basis;
53 using BackendType = typename DirichletValues::BackendType;
54 using FlagsType = typename DirichletValues::FlagsType;
55 using MultiIndex = typename Basis::MultiIndex;
56 using LocalView = typename Basis::LocalView;
57 using Intersection = typename Basis::GridView::Intersection;
58
59 pybind11::module scopedf = pybind11::module::import("dune.functions");
60 typedef Dune::Python::LocalViewWrapper<Basis> LocalViewWrapper;
61 auto includes = Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"};
62 auto lv = Dune::Python::insertClass<LocalViewWrapper>(
63 scopedf, "LocalView",
64 Dune::Python::GenerateTypeName("Dune::Python::LocalViewWrapper", Dune::MetaType<Basis>()), includes)
65 .first;
66
67 cls.def(pybind11::init([](const Basis& basis) { return new DirichletValues(basis); }),
68 pybind11::keep_alive<1, 2>());
69
70 // Eigen::Ref needed due to https://pybind11.readthedocs.io/en/stable/advanced/cast/eigen.html#pass-by-reference
71 cls.def("fixBoundaryDOFs",
72 [](DirichletValues& self, const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int)>& f) {
73 auto lambda = [&](BackendType& vec, const MultiIndex& indexGlobal) {
74 // we explicitly only allow flat indices
75 f(vec.vector(), indexGlobal[0]);
76 };
77 self.fixBoundaryDOFs(lambda);
78 });
79
80 cls.def("fixBoundaryDOFsUsingLocalView",
81 [](DirichletValues& self,
82 const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LocalViewWrapper&)>& f) {
83 auto lambda = [&](BackendType& vec, int localIndex, LocalView& lv) {
84 auto lvWrapper = LocalViewWrapper(lv.globalBasis());
85 // this can be simplified when
86 // https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418 becomes available
87 pybind11::object obj = pybind11::cast(lv.element());
88 lvWrapper.bind(obj);
89 f(vec.vector(), localIndex, lvWrapper);
90 };
91 self.fixBoundaryDOFs(lambda);
92 });
93
94 cls.def("fixBoundaryDOFsUsingLocalViewAndIntersection",
95 [](DirichletValues& self,
96 const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LocalViewWrapper&, const Intersection&)>&
97 f) {
98 auto lambda = [&](BackendType& vec, int localIndex, LocalView& lv, const Intersection& intersection) {
99 auto lvWrapper = LocalViewWrapper(lv.globalBasis());
100 // this can be simplified when
101 // https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418 becomes available
102 pybind11::object obj = pybind11::cast(lv.element());
103 lvWrapper.bind(obj);
104 f(vec.vector(), localIndex, lvWrapper, intersection);
105 };
106 self.fixBoundaryDOFs(lambda);
107 });
108
109 cls.def("fixDOFs",
110 [](DirichletValues& self, const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, const Basis&)>& f) {
111 self.fixBoundaryDOFs(f);
112 });
113 }
114
115} // 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:49
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:81
Definition: flatassembler.hh:20
def basis(gv, tree)
Definition: basis.py:10
Wrapper class for a hierarchical basis constructed from a pre-basis.
Definition: utils/basis.hh:29
Class for handling Dirichlet boundary conditions in Ikarus.
Definition: utils/dirichletvalues.hh:46
decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >())) BackendType
Definition: utils/dirichletvalues.hh:51
void fixBoundaryDOFs(F &&f)
Function to fix (set boolean values to true or false) degrees of freedom on the boundary.
Definition: utils/dirichletvalues.hh:67
FlagsType_ FlagsType
Definition: utils/dirichletvalues.hh:49
std::remove_cvref_t< Basis_ > Basis
Definition: utils/dirichletvalues.hh:48