version 0.4.1
registerpreelement.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#pragma once
5
6#include <dune/common/fvector.hh>
7#include <dune/python/pybind11/pybind11.h>
8
9#include <Eigen/Core>
10
11namespace Ikarus::Python {
12
22template <class NonLinearElasticPre, class... options>
23void registerNonLinearElasticPre(pybind11::handle scope, pybind11::class_<NonLinearElasticPre, options...> cls) {
25 cls.def(pybind11::init([](const Material& mat) { return new NonLinearElasticPre(mat); }));
26}
27
37template <class LinearElasticPre, class... options>
38void registerLinearElasticPre(pybind11::handle scope, pybind11::class_<LinearElasticPre, options...> cls) {
39 using Material = typename LinearElasticPre::Material;
40 cls.def(pybind11::init([](const Material& mat) { return new LinearElasticPre(mat); }));
41}
42
52template <class TrussPre, class... options>
53void registerTrussPre(pybind11::handle scope, pybind11::class_<TrussPre, options...> cls) {
54 cls.def(pybind11::init([](double emod, double area) { return new TrussPre({emod, area}); }));
55}
56
66template <class KirchhoffLoveShellPre, class... options>
67void registerKirchhoffLoveShellPre(pybind11::handle scope, pybind11::class_<KirchhoffLoveShellPre, options...> cls) {
68 cls.def(pybind11::init(
69 [](const double& E, const double& nu, const double& h) { return new KirchhoffLoveShellPre({E, nu}, h); }));
70}
71
81template <class EASPre, class... options>
82void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_<EASPre, options...> cls) {
83 cls.def(pybind11::init([](int numberOfParameter) { return new EASPre(numberOfParameter); }));
84}
85
95template <class NeumannBoundaryLoadPre, class... options>
96void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_<NeumannBoundaryLoadPre, options...> cls) {
97 using BoundaryPatchType = typename NeumannBoundaryLoadPre::BoundaryPatchType;
98 using GridView = typename NeumannBoundaryLoadPre::GridView;
99
100 using LoadFunction = std::function<Eigen::Vector<double, NeumannBoundaryLoadPre::worldDim>(
102 cls.def(pybind11::init([](const BoundaryPatchType& patch, LoadFunction volumeLoad) {
103 return new NeumannBoundaryLoadPre(&patch, volumeLoad);
104 }),
105 pybind11::keep_alive<1, 2>());
106}
107
117template <class VolumeLoadPre, class... options>
118void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_<VolumeLoadPre, options...> cls) {
119 using LoadFunction = std::function<Eigen::Vector<double, VolumeLoadPre::worldDim>(
121 cls.def(pybind11::init([](LoadFunction volumeLoad) { return new VolumeLoadPre(volumeLoad); }));
122}
123
124} // 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:127
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:38
void registerKirchhoffLoveShellPre(pybind11::handle scope, pybind11::class_< KirchhoffLoveShellPre, options... > cls)
Registers a KirchhoffLoveShellPre class in Python.
Definition: registerpreelement.hh:67
void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_< VolumeLoadPre, options... > cls)
Registers a VolumeLoadPre class in Python.
Definition: registerpreelement.hh:118
void registerTrussPre(pybind11::handle scope, pybind11::class_< TrussPre, options... > cls)
Registers a TrussPre class in Python.
Definition: registerpreelement.hh:53
void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_< EASPre, options... > cls)
Registers an EnhancedAssumedStrainsPre class in Python.
Definition: registerpreelement.hh:82
void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_< NeumannBoundaryLoadPre, options... > cls)
Registers a NeumannBoundaryLoadPre class in Python.
Definition: registerpreelement.hh:96
void registerNonLinearElasticPre(pybind11::handle scope, pybind11::class_< NonLinearElasticPre, options... > cls)
Registers a NonLinearElasticPre class in Python.
Definition: registerpreelement.hh:23
A PreFE struct for Kirchhoff-Love shell elements.
Definition: kirchhoffloveshell.hh:31
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:40
MAT Material
Definition: linearelastic.hh:41
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
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:43
MAT Material
Definition: nonlinearelastic.hh:44
A PreFE struct for truss elements.
Definition: truss.hh:28
Definition: utils/dirichletvalues.hh:32