version 0.3.7
python/finiteelements/nonlinearelastic.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
4#pragma once
5
6#include <dune/fufem/boundarypatch.hh>
7#include <dune/functions/functionspacebases/lagrangebasis.hh>
8#include <dune/functions/functionspacebases/powerbasis.hh>
9#include <dune/grid/yaspgrid.hh>
10#include <dune/python/common/typeregistry.hh>
11#include <dune/python/functions/globalbasis.hh>
12#include <dune/python/pybind11/eigen.h>
13#include <dune/python/pybind11/functional.h>
14#include <dune/python/pybind11/pybind11.h>
15#include <dune/python/pybind11/stl.h>
16
19#include <ikarus/utils/basis.hh>
20
21namespace Ikarus::Python {
22
23 template <class NonLinearElastic, class... options>
24 void registerNonLinearElastic(pybind11::handle scope, pybind11::class_<NonLinearElastic, options...> cls) {
25 using pybind11::operator""_a;
26
27 using GlobalBasis = typename NonLinearElastic::Basis;
28 using FlatBasis = typename NonLinearElastic::FlatBasis;
29 using GridView = typename GlobalBasis::GridView;
30 using Element = typename NonLinearElastic::Element;
31 using Traits = typename NonLinearElastic::Traits;
33 using Material = typename NonLinearElastic::Material;
35 using ResultRequirementsType = typename NonLinearElastic::ResultRequirementsType;
36
37 cls.def(pybind11::init([](const GlobalBasis& basis, const Element& element, const Material& mat) {
38 return new NonLinearElastic(basis, element, mat);
39 }),
40 pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>());
41
42 using LoadFunction = std::function<Eigen::Vector<double, Traits::worlddim>(
44
45 cls.def(pybind11::init(
46 [](const GlobalBasis& basis, const Element& element, const Material& mat,
47 const LoadFunction volumeLoad) { return new NonLinearElastic(basis, element, mat, volumeLoad); }),
48 pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>());
49
50 cls.def(pybind11::init([](const GlobalBasis& basis, const Element& element, const Material& mat,
51 const LoadFunction volumeLoad, const BoundaryPatch<GridView>& bp,
52 const LoadFunction neumannBoundaryLoad) {
53 return new NonLinearElastic(basis, element, mat, volumeLoad, &bp, neumannBoundaryLoad);
54 }),
55 pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>(), pybind11::keep_alive<1, 7>());
56
57 pybind11::module scopedf = pybind11::module::import("dune.functions");
58
59 typedef Dune::Python::LocalViewWrapper<FlatBasis> LocalViewWrapper;
60 auto includes = Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"};
61 auto lv = Dune::Python::insertClass<LocalViewWrapper>(
62 scopedf, "LocalViewWrapper",
63 Dune::Python::GenerateTypeName("Dune::Python::LocalViewWrapperWrapper", Dune::MetaType<FlatBasis>()),
64 includes)
65 .first;
66 lv.def("bind", &LocalViewWrapper::bind);
67 lv.def("unbind", &LocalViewWrapper::unbind);
68 lv.def("index", [](const LocalViewWrapper& localView, int index) { return localView.index(index); });
69 lv.def("__len__", [](LocalViewWrapper& self) -> int { return self.size(); });
70
71 Dune::Python::Functions::registerTree<typename LocalViewWrapper::Tree>(lv);
72 lv.def("tree", [](const LocalViewWrapper& view) { return view.tree(); });
73
74 auto basisName = Dune::className<FlatBasis>();
75 auto entry = Dune::Python::insertClass<FlatBasis>(
76 scopedf, "DefaultGlobalBasis", pybind11::buffer_protocol(), Dune::Python::GenerateTypeName(basisName),
77 Dune::Python::IncludeFiles{"dune/python/functions/globalbasis.hh"});
78 if (entry.second) Dune::Python::registerGlobalBasis(scopedf, entry.first);
79
80 cls.def(
81 "localView",
82 [](NonLinearElastic& self) -> LocalViewWrapper {
83 auto lvWrapped = LocalViewWrapper(self.localView().globalBasis());
84 // this can be simplified when https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/418
85 // becomes available
86 pybind11::object obj = pybind11::cast(self.localView().element());
87 lvWrapped.bind(obj);
88 return lvWrapped;
89 },
90 pybind11::keep_alive<0, 1>());
91 cls.def("calculateScalar",
92 [](NonLinearElastic& self, const FERequirements& req) { return self.calculateScalar(req); });
93 cls.def("calculateVector", [](NonLinearElastic& self, const FERequirements& req, Eigen::Ref<Eigen::VectorXd> vec) {
94 return self.calculateVector(req, vec);
95 });
96 cls.def(
97 "calculateMatrix",
98 [](NonLinearElastic& self, const FERequirements& req, Eigen::Ref<Eigen::MatrixXd> mat) {
99 return self.calculateMatrix(req, mat);
100 },
101 pybind11::arg("FERequirements"), pybind11::arg("elementMatrix").noconvert());
102
103 cls.def(
104 "resultAt",
105 [](NonLinearElastic& self, const ResultRequirementsType& req,
107 ResultTypeMap<double> resultTypeMap;
108 self.calculateAt(req, local, resultTypeMap);
109 if (resType == ResultType::noType)
110 return resultTypeMap.getSingleResult().second;
111 else
112 return resultTypeMap.getResult(resType);
113 },
114 pybind11::arg("resultRequirements"), pybind11::arg("local"), pybind11::arg("resultType") = ResultType::noType);
115 }
116
117} // namespace Ikarus::Python
Definition of the LinearElastic class for finite element mechanics computations.
ResultType
A strongly typed enum class representing the type of the result request.
Definition: ferequirements.hh:103
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:81
Definition: flatassembler.hh:20
void registerNonLinearElastic(pybind11::handle scope, pybind11::class_< NonLinearElastic, options... > cls)
Definition: python/finiteelements/nonlinearelastic.hh:24
def basis(gv, tree)
Definition: basis.py:10
def NonLinearElastic(basis, element, material, volumeLoad=None, bp=None, neumannBoundaryLoad=None)
Definition: finite_elements/__init__.py:149
Class representing the requirements for finite element calculations.
Definition: ferequirements.hh:170
Class representing a map of result types to result arrays.
Definition: ferequirements.hh:342
auto & getSingleResult()
Get the result array for a single result type.
Definition: ferequirements.hh:379
ResultArray & getResult(const ResultType &resultType)
Get the result array for a specific result type.
Definition: ferequirements.hh:368
Interface classf or materials.
Definition: interface.hh:75
NonLinearElastic class represents a non-linear elastic finite element.
Definition: finiteelements/mechanics/nonlinearelastic.hh:44
typename LocalView::Element Element
Definition: finiteelements/mechanics/nonlinearelastic.hh:53
TraitsFromLocalView< LocalView, useEigenRef > Traits
Definition: finiteelements/mechanics/nonlinearelastic.hh:56
Basis_ Basis
Definition: finiteelements/mechanics/nonlinearelastic.hh:46
FERequirements_ FERequirementType
Definition: finiteelements/mechanics/nonlinearelastic.hh:50
ResultRequirements< FERequirementType > ResultRequirementsType
Definition: finiteelements/mechanics/nonlinearelastic.hh:51
Material_ Material
Definition: finiteelements/mechanics/nonlinearelastic.hh:47
typename Basis::FlatBasis FlatBasis
Definition: finiteelements/mechanics/nonlinearelastic.hh:48
Definition: resultevaluators.hh:20
Definition of the NonLinearElastic class for finite element mechanics computations.
Wrapper around Dune-functions global basis.