version 0.4.1
registerpreelement.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/python/pybind11/pybind11.h>
7
8namespace Ikarus::Python {
9
19template <class NonLinearElasticPre, class... options>
20void registerNonLinearElasticPre(pybind11::handle scope, pybind11::class_<NonLinearElasticPre, options...> cls) {
22 cls.def(pybind11::init([](const Material& mat) { return new NonLinearElasticPre(mat); }));
23}
24
34template <class LinearElasticPre, class... options>
35void registerLinearElasticPre(pybind11::handle scope, pybind11::class_<LinearElasticPre, options...> cls) {
36 using Material = typename LinearElasticPre::Material;
37 cls.def(pybind11::init([](const Material& mat) { return new LinearElasticPre(mat); }));
38}
39
49template <class TrussPre, class... options>
50void registerTrussPre(pybind11::handle scope, pybind11::class_<TrussPre, options...> cls) {
51 cls.def(pybind11::init([](double emod, double area) { return new TrussPre({emod, area}); }));
52}
53
63template <class KirchhoffLoveShellPre, class... options>
64void registerKirchhoffLoveShellPre(pybind11::handle scope, pybind11::class_<KirchhoffLoveShellPre, options...> cls) {
65 cls.def(pybind11::init(
66 [](const double& E, const double& nu, const double& h) { return new KirchhoffLoveShellPre({E, nu}, h); }));
67}
68
78template <class EASPre, class... options>
79void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_<EASPre, options...> cls) {
80 cls.def(pybind11::init([](int numberOfParameter) { return new EASPre(numberOfParameter); }));
81}
82
92template <class NeumannBoundaryLoadPre, class... options>
93void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_<NeumannBoundaryLoadPre, options...> cls) {
94 using BoundaryPatchType = typename NeumannBoundaryLoadPre::BoundaryPatchType;
95 using GridView = typename NeumannBoundaryLoadPre::GridView;
96
97 using LoadFunction = std::function<Eigen::Vector<double, NeumannBoundaryLoadPre::worldDim>(
99 cls.def(pybind11::init([](const BoundaryPatchType& patch, LoadFunction volumeLoad) {
100 return new NeumannBoundaryLoadPre(&patch, volumeLoad);
101 }),
102 pybind11::keep_alive<1, 2>());
103}
104
114template <class VolumeLoadPre, class... options>
115void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_<VolumeLoadPre, options...> cls) {
116 using LoadFunction = std::function<Eigen::Vector<double, VolumeLoadPre::worldDim>(
118 cls.def(pybind11::init([](LoadFunction volumeLoad) { return new VolumeLoadPre(volumeLoad); }));
119}
120
121} // namespace Ikarus::Python
auto volumeLoad(const std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> &f)
A helper function to create a volume load skill.
Definition: volume.hh:128
VolumeLoadPre(F f) -> VolumeLoadPre< traits::FunctionTraits< F >::return_type::RowsAtCompileTime >
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:82
Definition: flatassembler.hh:21
void registerLinearElasticPre(pybind11::handle scope, pybind11::class_< LinearElasticPre, options... > cls)
Registers a LinearElasticPre class in Python.
Definition: registerpreelement.hh:35
void registerKirchhoffLoveShellPre(pybind11::handle scope, pybind11::class_< KirchhoffLoveShellPre, options... > cls)
Registers a KirchhoffLoveShellPre class in Python.
Definition: registerpreelement.hh:64
void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_< VolumeLoadPre, options... > cls)
Registers a VolumeLoadPre class in Python.
Definition: registerpreelement.hh:115
void registerTrussPre(pybind11::handle scope, pybind11::class_< TrussPre, options... > cls)
Registers a TrussPre class in Python.
Definition: registerpreelement.hh:50
void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_< EASPre, options... > cls)
Registers an EnhancedAssumedStrainsPre class in Python.
Definition: registerpreelement.hh:79
void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_< NeumannBoundaryLoadPre, options... > cls)
Registers a NeumannBoundaryLoadPre class in Python.
Definition: registerpreelement.hh:93
void registerNonLinearElasticPre(pybind11::handle scope, pybind11::class_< NonLinearElasticPre, options... > cls)
Registers a NonLinearElasticPre class in Python.
Definition: registerpreelement.hh:20
A PreFE struct for Kirchhoff-Love shell elements.
Definition: kirchhoffloveshell.hh:31
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:38
MAT Material
Definition: linearelastic.hh:39
A PreFE struct for Neumann boundary load skill.
Definition: traction.hh:23
GV GridView
Definition: traction.hh:24
BoundaryPatch< GridView > BoundaryPatchType
Definition: traction.hh:29
Interface classf or materials.
Definition: finiteelements/mechanics/materials/interface.hh:80
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:42
MAT Material
Definition: nonlinearelastic.hh:43
A PreFE struct for truss elements.
Definition: truss.hh:28
Definition: utils/dirichletvalues.hh:32