version 0.4.1
flatassembler.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/pybind11/eigen.h>
13#include <dune/python/pybind11/pybind11.h>
14#include <dune/python/pybind11/stl.h>
15#include <dune/python/pybind11/stl_bind.h>
16
18#include <ikarus/utils/basis.hh>
19
20namespace Ikarus::Python {
21
22#define MAKE_ASSEMBLER_REGISTERY_FUNCTION(name) \
23 \
44 template <class Assembler, class... options> \
45 void register##name(pybind11::handle scope, pybind11::class_<Assembler, options...> cls) { \
46 using pybind11::operator""_a; \
47 using FEContainer = typename Assembler::FEContainer; \
48 using Basis = typename Assembler::Basis; \
49 using DirichletValuesType = typename Assembler::DirichletValuesType; \
50 using FERequirementType = typename Assembler::FERequirementType; \
51 pybind11::module m = pybind11::module::import("ikarus"); \
52 cls.def(pybind11::init([](const pybind11::list& fes, const DirichletValuesType& dirichletValues) { \
53 /*here a copy of the whole vector of fes takes place! There is no way to prevent this if we want that \
54 * the user can pass native python lists here, see \
55 * https://pybind11.readthedocs.io/en/stable/advanced/cast/stl.html */ \
56 FEContainer fesV = fes.template cast<FEContainer>(); \
57 return new Assembler(std::move(fesV), dirichletValues); \
58 }), \
59 pybind11::keep_alive<1, 3>()); \
60 \
61 /* sparse matrices need to be copied to python therefore we remove the reference of the return type, see */ \
62 /* https://github.com/pybind/pybind11/blob/cbb876cc7b02c5f57e715cbc2c46ead3d1fbcf79/tests/test_eigen_matrix.cpp#L332-L341 \
63 */ \
64 cls.def( \
65 "getMatrix", \
66 [](Assembler& self, const FERequirementType& req) -> std::remove_cvref_t<decltype(self.getMatrix(req))> { \
67 return self.getMatrix(req); \
68 }, \
69 pybind11::return_value_policy::copy); \
70 cls.def( \
71 "getReducedMatrix", \
72 [](Assembler& self, const FERequirementType& req) \
73 -> std::remove_cvref_t<decltype(self.getReducedMatrix(req))> { return self.getReducedMatrix(req); }, \
74 pybind11::return_value_policy::copy); \
75 cls.def( \
76 "getVector", [](Assembler& self, const FERequirementType& req) { return self.getVector(req); }, \
77 pybind11::return_value_policy::reference); \
78 cls.def( \
79 "getScalar", [](Assembler& self, const FERequirementType& req) { return self.getScalar(req); }, \
80 pybind11::return_value_policy::copy); \
81 cls.def( \
82 "getReducedVector", [](Assembler& self, const FERequirementType& req) { return self.getReducedVector(req); }, \
83 pybind11::return_value_policy::reference); \
84 cls.def( \
85 "createFullVector", \
86 [](Assembler& self, Eigen::Ref<const Eigen::VectorXd> redVec) { return self.createFullVector(redVec); }, \
87 pybind11::return_value_policy::move); \
88 cls.def("reducedSize", [](Assembler& self) { return self.reducedSize(); }, pybind11::return_value_policy::copy); \
89 }
90
93
94} // namespace Ikarus::Python
#define MAKE_ASSEMBLER_REGISTERY_FUNCTION(name)
Definition: flatassembler.hh:22
Definition of the LinearElastic class for finite element mechanics computations.
Definition: flatassembler.hh:20
SparseFlatAssembler assembles matrix quantities using a flat basis Indexing strategy....
Definition: simpleassemblers.hh:277
DenseFlatAssembler assembles matrix quantities using a flat basis Indexing strategy....
Definition: simpleassemblers.hh:385
Wrapper around Dune-functions global basis.