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 | Basis = Basis_ |
| using | FlatBasis = typename Basis::FlatBasis |
| using | BasePowerFE = PowerBasisFE< FlatBasis > |
| using | FERequirementType = FERequirements_ |
| using | ResultRequirementsType = ResultRequirements< FERequirementType > |
| using | LocalView = typename FlatBasis::LocalView |
| using | Element = typename LocalView::Element |
| using | Geometry = typename Element::Geometry |
| using | GridView = typename FlatBasis::GridView |
| using | Traits = TraitsFromLocalView< LocalView, useEigenRef > |
| using | LocalBasisType = decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) |
| using | RootBasis = Basis_::FlatBasis |
| Type of the root basis. More... | |
| using | GlobalIndex = typename LocalView::MultiIndex |
| Type of the global index. More... | |
| using | GridElement = typename LocalView::Element |
| Type of the grid element. More... | |
Public Member Functions | |
| template<typename VolumeLoad = utils::LoadDefault, typename NeumannBoundaryLoad = utils::LoadDefault> | |
| KirchhoffLoveShell (const Basis &globalBasis, const typename LocalView::Element &element, double emod, double nu, double thickness, VolumeLoad p_volumeLoad={}, const BoundaryPatch< GridView > *p_neumannBoundary=nullptr, NeumannBoundaryLoad p_neumannBoundaryLoad={}) | |
| Constructor for the KirchhoffLoveShell class. More... | |
| template<typename ScalarType = double> | |
| auto | displacementFunction (const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const |
| Get the displacement function and nodal displacements. More... | |
| double | calculateScalar (const FERequirementType &req) const |
| Calculate the scalar value. More... | |
| void | calculateVector (const FERequirementType &req, typename Traits::template VectorType<> force) const |
| Calculate the vector associated with the given FERequirementType. More... | |
| void | calculateMatrix (const FERequirementType &req, typename Traits::template MatrixType<> K) const |
| Calculate the matrix associated with the given FERequirementType. More... | |
| void | calculateAt (const ResultRequirementsType &req, const Dune::FieldVector< double, Traits::mydim > &local, ResultTypeMap< double > &result) const |
| Calculate results at local coordinates. More... | |
| constexpr size_t | size () const |
| Get the size of the local view. More... | |
| void | globalFlatIndices (std::vector< GlobalIndex > &globalIndices) const |
| Get the global flat indices for the power basis. More... | |
| const GridElement & | gridElement () const |
| Get the grid element associated with the local view. More... | |
| const LocalView & | localView () const |
| Get the const reference to the local view. More... | |
| LocalView & | localView () |
| Get the reference to the local view. More... | |
Public Attributes | |
| std::shared_ptr< const Geometry > | geo_ |
| Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > | localBasis |
| std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> | volumeLoad |
| std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> | neumannBoundaryLoad |
| const BoundaryPatch< GridView > * | neumannBoundary |
| DefaultMembraneStrain | membraneStrain |
| double | emod_ |
| double | nu_ |
| double | thickness_ |
| size_t | numberOfNodes {0} |
| int | order {} |
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 |
| static constexpr int | num_children |
| Number of children in the powerBasis. More... | |
Protected Member Functions | |
| 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 ScalarType > | |
| void | calculateMatrixImpl (const FERequirementType &par, typename Traits::template MatrixType< ScalarType > K, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const |
| template<typename ScalarType > | |
| void | calculateVectorImpl (const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const |
| template<typename ScalarType > | |
| auto | calculateScalarImpl (const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType |
| template<typename ScalarType > | |
| Eigen::Matrix< ScalarType, 3, 3 > | kgBending (const Eigen::Matrix< ScalarType, 2, 3 > &jcur, const Eigen::Matrix3< ScalarType > &h, const auto &dN, const auto &ddN, const Eigen::Vector3< ScalarType > &a3N, const Eigen::Vector3< ScalarType > &a3, const Eigen::Vector3< ScalarType > &S, int I, int J) const |
| template<typename ScalarType > | |
| Eigen::Matrix< ScalarType, 3, 3 > | bopBending (const Eigen::Matrix< ScalarType, 2, 3 > &jcur, const Eigen::Matrix3< ScalarType > &h, const auto &dN, const auto &ddN, const int node, const Eigen::Vector3< ScalarType > &a3N, const Eigen::Vector3< ScalarType > &a3) const |
| Eigen::Matrix< double, 3, 3 > | materialTangent (const Eigen::Matrix< double, 3, 3 > &Aconv) const |
This class represents Kirchhoff-Love shell finite elements.
| Basis_ | The basis type for the finite element. |
| FERequirements_ | The type representing the requirements for finite element calculations. |
| useEigenRef | A boolean indicating whether to use Eigen references for efficiency. |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::BasePowerFE = PowerBasisFE<FlatBasis> |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::Basis = Basis_ |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::Element = typename LocalView::Element |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::FERequirementType = FERequirements_ |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::FlatBasis = typename Basis::FlatBasis |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::Geometry = typename Element::Geometry |
|
inherited |
|
inherited |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::GridView = typename FlatBasis::GridView |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis()) |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::LocalView = typename FlatBasis::LocalView |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::ResultRequirementsType = ResultRequirements<FERequirementType> |
|
inherited |
| using Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::Traits = TraitsFromLocalView<LocalView, useEigenRef> |
|
inline |
Initializes the KirchhoffLoveShell instance with the given parameters.
| VolumeLoad | The type representing the volume load function. |
| NeumannBoundaryLoad | The type representing the Neumann boundary load function. |
| globalBasis | The global basis for the finite element. |
| element | The local element to bind. |
| emod | Young's modulus of the material. |
| nu | Poisson's ratio of the material. |
| thickness | Thickness of the shell. |
| p_volumeLoad | The volume load function (optional, default is utils::LoadDefault). |
| p_neumannBoundary | The Neumann boundary patch (optional, default is nullptr). |
| p_neumannBoundaryLoad | The Neumann boundary load function (optional, default is LoadDefault). |
|
inlineprotected |
|
inline |
Calculates the results at the specified local coordinates based on the given requirements and stores them in the result container.
| req | The result requirements. |
| local | The local coordinates at which results are to be calculated. |
| result | The result container to store the calculated values. |
|
inline |
| ScalarType | The scalar type for the calculation. |
| req | The FERequirementType object specifying the requirements for the calculation. |
| K | The matrix to store the calculated result. |
|
inlineprotected |
|
inline |
Calculates the scalar value based on the given FERequirements.
| req | The FERequirements. |
|
inlineprotected |
|
inline |
| ScalarType | The scalar type for the calculation. |
| req | The FERequirementType object specifying the requirements for the calculation. |
| force | The vector to store the calculated result. |
|
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. |
| gpPos | The type of the integration point position. |
| gpIndex | The type of the integration point index. |
| geo | The type of the geometry object. |
| uFunction | The type of the displacement field function. |
|
inline |
Retrieves the displacement function and nodal displacements based on the given FERequirements.
| ScalarType | The scalar type used for calculations. |
| par | The FERequirements. |
| dx | Optional additional displacement vector. |
|
inlineinherited |
The global indices are collected in a FlatInterLeaved order. I.e. x_0, y_0, z_0, ..., x_n, y_n, z_n
| globalIndices | Output vector to store global indices. |
|
inlineinherited |
|
inlineprotected |
|
inlineinherited |
|
inlineinherited |
|
inlineprotected |
|
inlineconstexprinherited |
|
staticconstexpr |
| double Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::emod_ |
| std::shared_ptr<const Geometry> Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::geo_ |
| Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType> > Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::localBasis |
| DefaultMembraneStrain Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::membraneStrain |
|
staticconstexpr |
|
staticconstexpr |
| const BoundaryPatch<GridView>* Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::neumannBoundary |
| std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)> Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::neumannBoundaryLoad |
| double Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::nu_ |
|
staticconstexprinherited |
| size_t Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::numberOfNodes {0} |
| int Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::order {} |
| double Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::thickness_ |
| std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)> Ikarus::KirchhoffLoveShell< Basis_, FERequirements_, useEigenRef >::volumeLoad |
|
staticconstexpr |