version 0.4.1
material.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/common/classname.hh>
12#include <dune/python/common/typeregistry.hh>
13#include <dune/python/pybind11/eigen.h>
14#include <dune/python/pybind11/pybind11.h>
15#include <dune/python/pybind11/stl.h>
16
19
20#define MAKE_MaterialFunction(clsName, materialName, functionname, vecSize) \
21 clsName.def( \
22 #functionname, \
23 [](materialName& self, const std::string& straintag, Eigen::Ref<const Eigen::Vector<double, vecSize>> eVoigt_) { \
24 if constexpr (not Concepts::IsMaterial<LinearElasticityT, materialName>) { \
25 Eigen::Vector<double, vecSize> eVoigt = eVoigt_; \
26 if (straintag == toString(StrainTags::rightCauchyGreenTensor)) \
27 return self.template functionname<StrainTags::rightCauchyGreenTensor>(eVoigt); \
28 else if (straintag == toString(StrainTags::greenLagrangian)) \
29 return self.template functionname<StrainTags::greenLagrangian>(eVoigt); \
30 else if (straintag == toString(StrainTags::linear)) \
31 DUNE_THROW(Dune::MathError, "Passing linear strain to " + std::string(#materialName) + \
32 " does not makes sense use LinearElastic class"); \
33 else if (straintag == toString(StrainTags::displacementGradient)) \
34 DUNE_THROW(Dune::MathError, \
35 "Passing displacementGradient strain in 6d Voigt notation does not make any sense!"); \
36 else if (straintag == toString(StrainTags::deformationGradient)) \
37 DUNE_THROW(Dune::MathError, \
38 "Passing deformationGradient strain in 6d Voigt notation does not make any sense!"); \
39 else \
40 DUNE_THROW(Dune::MathError, straintag + "is not a valid strain tag."); \
41 } else { \
42 Eigen::Vector<double, vecSize> eVoigt = eVoigt_; /* Linear elastic path */ \
43 if (straintag == toString(StrainTags::linear)) \
44 return self.template functionname<StrainTags::linear>(eVoigt); \
45 else \
46 DUNE_THROW(Dune::MathError, "Linear elastic material only accepts linear strains!"); \
47 } \
48 __builtin_unreachable(); \
49 }, \
50 "StrainName"_a, "strainVector"_a);
51
52#define MAKE_MATERIAL_REGISTERY_FUNCTION(Materialname, vecSize) \
53 template <class Materialname, class... options> \
54 void register##Materialname(pybind11::handle scope, pybind11::class_<Materialname, options...> cls##Materialname) { \
55 using pybind11::operator""_a; \
56 namespace py = pybind11; \
57 cls##Materialname.def(pybind11::init([](double emod, double nu) { \
58 auto matParameter = \
59 Ikarus::toLamesFirstParameterAndShearModulus({.emodul = emod, .nu = nu}); \
60 return new Materialname(matParameter); \
61 }), \
62 "Material constructor that takes Young's modulus E and Poisson's ratio nu", "E"_a, "nu"_a); \
63 MAKE_MaterialFunction(cls##Materialname, Materialname, storedEnergy, vecSize); \
64 MAKE_MaterialFunction(cls##Materialname, Materialname, stresses, vecSize); \
65 MAKE_MaterialFunction(cls##Materialname, Materialname, tangentModuli, vecSize); \
66 \
67 using PlaneStressClass = decltype(planeStress(std::declval<Materialname>())); \
68 auto includes = Dune::Python::IncludeFiles{"ikarus/finiteelements/mechanics/materials.hh"}; \
69 auto pS = Dune::Python::insertClass<PlaneStressClass>( \
70 scope, std::string("PlaneStress_") + #Materialname, \
71 Dune::Python::GenerateTypeName( \
72 "Ikarus::VanishingStress<std::array<Ikarus::Impl::StressIndexPair, " \
73 "3ul>{{Ikarus::Impl::StressIndexPair{2ul, 1ul}, Ikarus::Impl::StressIndexPair{2ul,0ul}, " \
74 "Ikarus::Impl::StressIndexPair{2ul, 2ul}}}," + \
75 Dune::className<Materialname>() + ">"), \
76 includes) \
77 .first; \
78 MAKE_MaterialFunction(pS, PlaneStressClass, storedEnergy, 3); \
79 MAKE_MaterialFunction(pS, PlaneStressClass, stresses, 3); \
80 MAKE_MaterialFunction(pS, PlaneStressClass, tangentModuli, 3); \
81 MAKE_MaterialFunction(pS, PlaneStressClass, storedEnergy, 6); \
82 MAKE_MaterialFunction(pS, PlaneStressClass, stresses, 6); \
83 MAKE_MaterialFunction(pS, PlaneStressClass, tangentModuli, 6); \
84 cls##Materialname.def("asPlaneStress", [](Materialname& self) { \
85 return planeStress(self); \
86 }); /* no keep_alive since planeStress copies the material */ \
87 using shellMaterialClass = decltype(shellMaterial(std::declval<Materialname>())); \
88 auto shellmaterial = \
89 Dune::Python::insertClass<shellMaterialClass>( \
90 scope, std::string("Shell_") + #Materialname, \
91 Dune::Python::GenerateTypeName("Ikarus::VanishingStress<std::array<Ikarus::Impl::StressIndexPair, " \
92 "1ul>{{Ikarus::Impl::StressIndexPair{2ul, 2ul}}}," + \
93 Dune::className<Materialname>() + ">"), \
94 includes) \
95 .first; \
96 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, storedEnergy, 5); \
97 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, stresses, 5); \
98 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, tangentModuli, 5); \
99 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, storedEnergy, 6); \
100 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, stresses, 6); \
101 MAKE_MaterialFunction(shellmaterial, shellMaterialClass, tangentModuli, 6); \
102 cls##Materialname.def("asShellMaterial", [](Materialname& self) { \
103 return shellMaterial(self); \
104 }); /* no keep_alive since shellMaterial copies the material */ \
105 using beamMaterialClass = decltype(beamMaterial(std::declval<Materialname>())); \
106 auto beammaterial = Dune::Python::insertClass<beamMaterialClass>( \
107 scope, std::string("Beam_") + #Materialname, \
108 Dune::Python::GenerateTypeName( \
109 "Ikarus::VanishingStress<std::array<Ikarus::Impl::StressIndexPair, " \
110 "2ul>{{Impl::StressIndexPair{1, 1},Ikarus::Impl::StressIndexPair{2ul, 2ul}}}," + \
111 Dune::className<Materialname>() + ">"), \
112 includes) \
113 .first; \
114 MAKE_MaterialFunction(beammaterial, beamMaterialClass, storedEnergy, 4); \
115 MAKE_MaterialFunction(beammaterial, beamMaterialClass, stresses, 4); \
116 MAKE_MaterialFunction(beammaterial, beamMaterialClass, tangentModuli, 4); \
117 MAKE_MaterialFunction(beammaterial, beamMaterialClass, storedEnergy, 6); \
118 MAKE_MaterialFunction(beammaterial, beamMaterialClass, stresses, 6); \
119 MAKE_MaterialFunction(beammaterial, beamMaterialClass, tangentModuli, 6); \
120 cls##Materialname.def("asBeamMaterial", [](Materialname& self) { \
121 return beamMaterial(self); \
122 }); /* no keep_alive since beamMaterial copies the material */ \
123 }
124
125#define MAKE_MATERIAL_CLASS_IN_MODULE(Materialname, args) \
126 auto includes##Materialname = Dune::Python::IncludeFiles{"ikarus/finiteelements/mechanics/materials.hh"}; \
127 auto cls##Materialname = \
128 Dune::Python::insertClass<Ikarus::Materialname<args>>( \
129 m, #Materialname, Dune::Python::GenerateTypeName("Ikarus::" + std::string(#Materialname<##args>)), \
130 includes##Materialname) \
131 .first; \
132 cls##Materialname.def(pybind11::init([](double emod, double nu) { \
133 auto matParameter = \
134 Ikarus::toLamesFirstParameterAndShearModulus({.emodul = emod, .nu = nu}); \
135 return new Materialname(matParameter); \
136 }), \
137 "Material constructor that takes Young's modulus E and Poisson's ratio nu" \
138 "E"_a, \
139 "nu"_a); \
140 MAKE_MaterialFunction(Materialname<##args>, storedEnergy); \
141 MAKE_MaterialFunction(Materialname<##args>, stresses); \
142 MAKE_MaterialFunction(Materialname<##args>, tangentModuli);
143
144namespace Ikarus::Python {
145
149} // namespace Ikarus::Python
Several concepts.
Header file for material models in Ikarus finite element mechanics.
Definition: flatassembler.hh:20
MAKE_MATERIAL_REGISTERY_FUNCTION(LinearElasticity, 6)
Implementation of the Linear Elasticity material model.The energy is computed as.
Definition: linearelasticity.hh:36
Implementation of the Neo-Hookean material model.The energy is computed as.
Definition: neohooke.hh:37
Implementation of the Saint Venant-Kirchhoff material model.The energy is computed as.
Definition: svk.hh:37