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 cls.def(pybind11::init([](double emod, double nu) { return new LinearElasticPre({emod, nu}); }));
37}
38
48template <class KirchhoffLoveShellPre, class... options>
49void registerKirchhoffLoveShellPre(pybind11::handle scope, pybind11::class_<KirchhoffLoveShellPre, options...> cls) {
50 cls.def(pybind11::init(
51 [](const double& E, const double& nu, const double& h) { return new KirchhoffLoveShellPre({E, nu}, h); }));
52}
53
63template <class EASPre, class... options>
64void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_<EASPre, options...> cls) {
65 cls.def(pybind11::init([](int numberOfParameter) { return new EASPre(numberOfParameter); }));
66}
67
77template <class NeumannBoundaryLoadPre, class... options>
78void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_<NeumannBoundaryLoadPre, options...> cls) {
79 using BoundaryPatchType = typename NeumannBoundaryLoadPre::BoundaryPatchType;
80 using GridView = typename NeumannBoundaryLoadPre::GridView;
81
82 using LoadFunction = std::function<Eigen::Vector<double, NeumannBoundaryLoadPre::worldDim>(
84 cls.def(pybind11::init([](const BoundaryPatchType& patch, LoadFunction volumeLoad) {
85 return new NeumannBoundaryLoadPre(&patch, volumeLoad);
86 }),
87 pybind11::keep_alive<1, 2>());
88}
89
99template <class VolumeLoadPre, class... options>
100void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_<VolumeLoadPre, options...> cls) {
101 using LoadFunction = std::function<Eigen::Vector<double, VolumeLoadPre::worldDim>(
103 cls.def(pybind11::init([](LoadFunction volumeLoad) { return new VolumeLoadPre(volumeLoad); }));
104}
105
106} // 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:20
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:49
void registerVolumeLoadPre(pybind11::handle scope, pybind11::class_< VolumeLoadPre, options... > cls)
Registers a VolumeLoadPre class in Python.
Definition: registerpreelement.hh:100
void registerEnhancedAssumedStrainsPre(pybind11::handle scope, pybind11::class_< EASPre, options... > cls)
Registers an EnhancedAssumedStrainsPre class in Python.
Definition: registerpreelement.hh:64
void registerNeumannBoundaryLoadPre(pybind11::handle scope, pybind11::class_< NeumannBoundaryLoadPre, options... > cls)
Registers a NeumannBoundaryLoadPre class in Python.
Definition: registerpreelement.hh:78
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:30
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:35
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: interface.hh:77
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:42
MAT Material
Definition: nonlinearelastic.hh:43
Definition: utils/dirichletvalues.hh:30