Kirchhoff-Love shell finite element class. More...
#include <ikarus/finiteelements/mechanics/kirchhoffloveshell.hh>
Classes | |
struct | KinematicVariables |
A structure representing kinematic variables. More... | |
Public Types | |
using | Traits = PreFE::Traits |
using | BasisHandler = typename Traits::BasisHandler |
using | FlatBasis = typename Traits::FlatBasis |
using | Requirement = FERequirementsFactory< FESolutions::displacement, FEParameter::loadfactor, Traits::useEigenRef >::type |
using | LocalView = typename Traits::LocalView |
using | Geometry = typename Traits::Geometry |
using | GridView = typename Traits::GridView |
using | Element = typename Traits::Element |
using | LocalBasisType = decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) |
using | Pre = KirchhoffLoveShellPre |
using | SupportedResultTypes = std::tuple< decltype(makeRT< ResultTypes >())... > |
Public Member Functions | |
KirchhoffLoveShell (const Pre &pre) | |
Constructor for the KirchhoffLoveShell class. More... | |
template<typename ST = double> | |
auto | displacementFunction (const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const |
Get the displacement function and nodal displacements. More... | |
Geometry | geometry () const |
size_t | numberOfNodes () const |
int | order () const |
template<template< typename, int, int > class RT> requires (supportsResultType<RT>()) | |
auto | calculateAtImpl (const Requirement &req, const Dune::FieldVector< double, Traits::mydim > &local) -> ResultWrapper< RT< double, myDim, worldDim >, ResultShape::Vector > |
Calculates a requested result at a specific local position. More... | |
Static Public Member Functions | |
static consteval bool | supportsResultType () |
Returns whether a ResultType is provided by the element. More... | |
Static Public Attributes | |
static constexpr int | myDim = Traits::mydim |
static constexpr int | worldDim = Traits::worlddim |
static constexpr int | membraneStrainSize = 3 |
static constexpr int | bendingStrainSize = 3 |
Protected Member Functions | |
void | bindImpl () |
A helper function to bind the local view to the element. More... | |
auto | computeMaterialAndStrains (const Dune::FieldVector< double, 2 > &gpPos, int gpIndex, const Geometry &geo, const auto &uFunction) const |
Compute material properties and strains at a given integration point. More... | |
template<typename ST > | |
void | calculateMatrixImpl (const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType< ST > K, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const |
template<typename ST > | |
void | calculateVectorImpl (const Requirement &par, const VectorAffordance &affordance, typename Traits::template VectorType< ST > force, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const |
template<typename ST > | |
auto | calculateScalarImpl (const Requirement &par, const ScalarAffordance &affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const -> ST |
template<typename ST > | |
Eigen::Matrix< ST, 3, 3 > | kgBending (const Eigen::Matrix< ST, 2, 3 > &jcur, const Eigen::Matrix3< ST > &h, const auto &dN, const auto &ddN, const Eigen::Vector3< ST > &a3N, const Eigen::Vector3< ST > &a3, const Eigen::Vector3< ST > &S, int I, int J) const |
template<typename ST > | |
Eigen::Matrix< ST, 3, 3 > | bopBending (const Eigen::Matrix< ST, 2, 3 > &jcur, const Eigen::Matrix3< ST > &h, const auto &dN, const auto &ddN, const int node, const Eigen::Vector3< ST > &a3N, const Eigen::Vector3< ST > &a3) const |
Eigen::Matrix< double, 3, 3 > | materialTangent (const Eigen::Matrix< double, 3, 3 > &Aconv) const |
Gets the material tangent matrix for the linear elastic material. More... | |
This class represents Kirchhoff-Love shell finite elements.
using Ikarus::KirchhoffLoveShell< PreFE, FE >::BasisHandler = typename Traits::BasisHandler |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::Element = typename Traits::Element |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::FlatBasis = typename Traits::FlatBasis |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::Geometry = typename Traits::Geometry |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::GridView = typename Traits::GridView |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis()) |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::LocalView = typename Traits::LocalView |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::Pre = KirchhoffLoveShellPre |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::Requirement = FERequirementsFactory<FESolutions::displacement, FEParameter::loadfactor, Traits::useEigenRef>::type |
|
inherited |
using Ikarus::KirchhoffLoveShell< PreFE, FE >::Traits = PreFE::Traits |
|
inline |
Initializes the KirchhoffLoveShell instance with the given parameters.
pre | The pre fe |
|
inlineprotected |
|
inlineprotected |
|
inline |
req | The FERequirementType object holding the global solution. |
local | Local position vector. |
RT | The type representing the requested result. |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
gpPos | The position of the integration point. |
gpIndex | The index of the integration point. |
geo | The geometry object providing position and derivatives. |
uFunction | The function representing the displacement field. |
|
inline |
Retrieves the displacement function and nodal displacements based on the given FERequirements.
ST | The scalar type used for calculations. |
par | The FERequirements. |
dx | Optional additional displacement vector. |
|
inline |
|
inlineprotected |
|
inlineprotected |
Aconv | A tensor to get the material tangent in the curvilinear coordinate system. |
|
inline |
|
inline |
|
inlinestaticinherited |
RT | requested ResultType |
|
staticconstexpr |
|
staticconstexpr |
|
staticconstexpr |
|
staticconstexpr |