4#include "../pythonhelpers.hh"
6#include <dune/python/pybind11/eigen.h>
7#include <dune/python/pybind11/operators.h>
8#include <dune/python/pybind11/pybind11.h>
10#include <spdlog/spdlog.h>
17 namespace py = pybind11;
18 using namespace pybind11::literals;
20 using namespace Eigen;
22 auto materials = m.def_submodule(
"materials",
"This is the submodule for materials in Ikarus");
28 pybind11::class_<Materials::LinearElasticity> linElastic(materials,
"LinearElasticity");
29 Ikarus::Python::registerLinearElasticity(materials, linElastic);
31 pybind11::class_<Materials::StVenantKirchhoff> svk(materials,
"StVenantKirchhoff");
32 Ikarus::Python::registerStVenantKirchhoff(materials, svk);
34 pybind11::class_<Materials::NeoHooke> nh(materials,
"NeoHooke");
35 Ikarus::Python::registerNeoHooke(materials, nh);
50 auto callTransformStrain =
52 if (eRaw.cols() == 1) {
53 Eigen::Vector<double, 6> E = eRaw;
54 return transformStrain<from_, to_>(E);
56 Eigen::Matrix<double, 3, 3> E = eRaw;
57 return transformStrain<from_, to_>(E);
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");
64 if (from == StrainTags::linear) {
65 spdlog::warn(
"No useful transformation available for linear strains");
71 if (to == StrainTags::greenLagrangian) {
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);
83 __builtin_unreachable();
85 }
else if (to == StrainTags::deformationGradient) {
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);
97 __builtin_unreachable();
99 }
else if (to == StrainTags::rightCauchyGreenTensor) {
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);
111 __builtin_unreachable();
114 __builtin_unreachable();
116 py::arg(
"to"), py::arg(
"from"), py::arg(
"strain"));
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.