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 <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
48template <class DirichletValues, class... options>
49void 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); }), pybind11::keep_alive<1, 2>());
68
69 // Eigen::Ref needed due to https://pybind11.readthedocs.io/en/stable/advanced/cast/eigen.html#pass-by-reference
70 cls.def("fixBoundaryDOFs",
71 [](DirichletValues& self, const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int)>& f) {
72 auto lambda = [&](BackendType& vec, const MultiIndex& indexGlobal) {
73 // we explicitly only allow flat indices
74 f(vec.vector(), indexGlobal[0]);
75 };
76 self.fixBoundaryDOFs(lambda);
77 });
78
79 cls.def("fixBoundaryDOFsUsingLocalView",
80 [](DirichletValues& self,
81 const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LocalViewWrapper&)>& f) {
82 auto lambda = [&](BackendType& vec, int localIndex, LocalView& lv) {
83 auto lvWrapper = LocalViewWrapper(lv.globalBasis());
84 // this can be simplified when
85 // https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418 becomes available
86 pybind11::object obj = pybind11::cast(lv.element());
87 lvWrapper.bind(obj);
88 f(vec.vector(), localIndex, lvWrapper);
89 };
90 self.fixBoundaryDOFs(lambda);
91 });
92
93 cls.def(
94 "fixBoundaryDOFsUsingLocalViewAndIntersection",
95 [](DirichletValues& self,
96 const std::function<void(Eigen::Ref<Eigen::VectorX<bool>>, int, LocalViewWrapper&, const Intersection&)>& f) {
97 auto lambda = [&](BackendType& vec, int localIndex, LocalView& lv, const Intersection& intersection) {
98 auto lvWrapper = LocalViewWrapper(lv.globalBasis());
99 // this can be simplified when
100 // https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418 becomes available
101 pybind11::object obj = pybind11::cast(lv.element());
102 lvWrapper.bind(obj);
103 f(vec.vector(), localIndex, lvWrapper, intersection);
104 };
105 self.fixBoundaryDOFs(lambda);
106 });
107
108 cls.def("fixDOFs",
109 [](DirichletValues& self, const std::function<void(const Basis&, Eigen::Ref<Eigen::VectorX<bool>>)>& f) {
110 auto lambda = [&](const Basis& basis, BackendType& vec) {
111 // we explicitly only allow flat indices
112 f(basis, vec.vector());
113 };
114 self.fixDOFs(lambda);
115 });
116}
117
118} // 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:82
Definition: flatassembler.hh:20
def basis(gv, tree)
Definition: basis.py:9
Class for handling Dirichlet boundary conditions in Ikarus.
Definition: utils/dirichletvalues.hh:47
std::remove_cvref_t< B > Basis
Definition: utils/dirichletvalues.hh:49
decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >())) BackendType
Definition: utils/dirichletvalues.hh:52
void fixBoundaryDOFs(F &&f)
Function to fix (set boolean values to true or false) degrees of freedom on the boundary.
Definition: utils/dirichletvalues.hh:70
void fixDOFs(F &&f)
Function to fix (set boolean values to true or false) degrees of freedom.
Definition: utils/dirichletvalues.hh:95
FC FlagsType
Definition: utils/dirichletvalues.hh:50