version 0.4.1
python/ikarus/materials/materials.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#include "../pythonhelpers.hh"
5
6#include <dune/python/pybind11/eigen.h>
7#include <dune/python/pybind11/operators.h>
8#include <dune/python/pybind11/pybind11.h>
9
10#include <spdlog/spdlog.h>
11
15
16void addMaterialsSubModule(pybind11::module& m) {
17 namespace py = pybind11;
18 using namespace pybind11::literals;
19 using namespace Ikarus;
20 using namespace Eigen;
21
22 auto materials = m.def_submodule("materials", "This is the submodule for materials in Ikarus");
23
27
28 pybind11::class_<Materials::LinearElasticity> linElastic(materials, "LinearElasticity");
29 Ikarus::Python::registerLinearElasticity(materials, linElastic);
30
31 pybind11::class_<Materials::StVenantKirchhoff> svk(materials, "StVenantKirchhoff");
32 Ikarus::Python::registerStVenantKirchhoff(materials, svk);
33
34 pybind11::class_<Materials::NeoHooke> nh(materials, "NeoHooke");
35 Ikarus::Python::registerNeoHooke(materials, nh);
36
47 materials.def(
48 "transformStrain",
49 [](StrainTags from, StrainTags to, Eigen::MatrixXd E) -> Eigen::MatrixXd {
50 auto callTransformStrain =
51 []<StrainTags from_, StrainTags to_>(const Eigen::MatrixXd& eRaw) -> Eigen::MatrixXd {
52 if (eRaw.cols() == 1) {
53 Eigen::Vector<double, 6> E = eRaw;
54 return transformStrain<from_, to_>(E);
55 } else {
56 Eigen::Matrix<double, 3, 3> E = eRaw;
57 return transformStrain<from_, to_>(E);
58 }
59 };
60 if (!((E.rows() == 3 && E.cols() == 3) || (E.rows() == 6 && E.cols() == 1)))
61 DUNE_THROW(Dune::IOError,
62 "Strain conversions are only implemented for matrices of dimension 3 or Voigt vectors of size 6");
63
64 if (from == StrainTags::linear) {
65 spdlog::warn("No useful transformation available for linear strains");
66 return E;
67 }
68 if (from == to)
69 return E;
70
71 if (to == StrainTags::greenLagrangian) {
72 switch (from) {
73 case StrainTags::deformationGradient:
74 return callTransformStrain
75 .template operator()<StrainTags::deformationGradient, StrainTags::greenLagrangian>(E);
77 return callTransformStrain
78 .template operator()<StrainTags::displacementGradient, StrainTags::greenLagrangian>(E);
80 return callTransformStrain
81 .template operator()<StrainTags::rightCauchyGreenTensor, StrainTags::greenLagrangian>(E);
82 default:
83 __builtin_unreachable();
84 }
85 } else if (to == StrainTags::deformationGradient) {
86 switch (from) {
87 case StrainTags::greenLagrangian:
88 return callTransformStrain
89 .template operator()<StrainTags::greenLagrangian, StrainTags::deformationGradient>(E);
91 return callTransformStrain
92 .template operator()<StrainTags::displacementGradient, StrainTags::deformationGradient>(E);
94 return callTransformStrain
95 .template operator()<StrainTags::rightCauchyGreenTensor, StrainTags::deformationGradient>(E);
96 default:
97 __builtin_unreachable();
98 }
99 } else if (to == StrainTags::rightCauchyGreenTensor) {
100 switch (from) {
101 case StrainTags::greenLagrangian:
102 return callTransformStrain
103 .template operator()<StrainTags::greenLagrangian, StrainTags::rightCauchyGreenTensor>(E);
105 return callTransformStrain
106 .template operator()<StrainTags::displacementGradient, StrainTags::rightCauchyGreenTensor>(E);
108 return callTransformStrain
109 .template operator()<StrainTags::deformationGradient, StrainTags::rightCauchyGreenTensor>(E);
110 default:
111 __builtin_unreachable();
112 }
113 }
114 __builtin_unreachable();
115 },
116 py::arg("to"), py::arg("from"), py::arg("strain"));
117}
Python bindings for materials.
Definition of several material related enums.
#define ENUM_BINDINGS_WITH_MODULE(Type, module)
Definition: pythonhelpers.hh:6
void addMaterialsSubModule(pybind11::module &m)
Definition: python/ikarus/materials/materials.hh:16
TangentModuliTags
A strongly typed enum class representing the type of the computed tangent moduli.
Definition: tags.hh:33
StressTags
A strongly typed enum class representing the type of the computed stresses.
Definition: tags.hh:26
StrainTags
A strongly typed enum class representing the type of the passed strain.
Definition: tags.hh:19
Definition: assemblermanipulatorbuildingblocks.hh:22
Definition: truncatedconjugategradient.hh:24
Header file for material models in Ikarus finite element mechanics.