LinearElastic class represents a linear elastic finite element. More...
#include <ikarus/finiteelements/mechanics/linearelastic.hh>
Public Types | |
using | Basis = Basis_ |
using | FlatBasis = typename Basis::FlatBasis |
using | BaseDisp = PowerBasisFE< FlatBasis > |
using | FERequirementType = FERequirements_ |
using | ResultRequirementsType = ResultRequirements< FERequirementType > |
using | LocalView = typename FlatBasis::LocalView |
using | Element = typename LocalView::Element |
using | Traits = TraitsFromLocalView< LocalView, useEigenRef > |
using | Geometry = typename Element::Geometry |
using | GridView = typename FlatBasis::GridView |
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> | |
LinearElastic (const Basis &globalBasis, const typename LocalView::Element &element, double emod, double nu, VolumeLoad p_volumeLoad={}, const BoundaryPatch< GridView > *p_neumannBoundary=nullptr, NeumannBoundaryLoad p_neumannBoundaryLoad={}) | |
Constructor for the LinearElastic class. More... | |
template<typename ScalarType = double> | |
auto | displacementFunction (const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const |
Gets the displacement function for the given FERequirementType and optional displacement vector. More... | |
template<class ScalarType = double> | |
auto | strainFunction (const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const |
Gets the strain function for the given FERequirementType and optional displacement vector. More... | |
auto | materialTangent () const |
Gets the material tangent matrix for the linear elastic material. More... | |
auto | materialTangentFunction (const FERequirementType &par) const |
Gets the material tangent function for the given FERequirementType. More... | |
double | calculateScalar (const FERequirementType &par) const |
Calculates the scalar energy for the given FERequirementType. More... | |
void | calculateVector (const FERequirementType &par, typename Traits::template VectorType<> force) const |
Calculates the vector force for the given FERequirementType. More... | |
void | calculateMatrix (const FERequirementType &par, typename Traits::template MatrixType<> K) const |
Calculates the matrix stiffness for the given FERequirementType. More... | |
void | calculateAt (const ResultRequirementsType &req, const Dune::FieldVector< double, Traits::mydim > &local, ResultTypeMap< double > &result) const |
Calculates results at a specific local position. 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 |
double | emod_ |
double | nu_ |
size_t | numberOfNodes {0} |
int | order {} |
Static Public Attributes | |
static constexpr int | myDim = Traits::mydim |
static constexpr int | worldDim = Traits::worlddim |
static constexpr int | num_children |
Number of children in the powerBasis. More... | |
Protected Member Functions | |
template<typename ScalarType > | |
auto | calculateScalarImpl (const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType |
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 |
Basis_ | The basis type for the finite element. |
FERequirements_ | The requirements for the finite element. |
useEigenRef | A boolean flag indicating whether to use Eigen references. |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::BaseDisp = PowerBasisFE<FlatBasis> |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::Basis = Basis_ |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::Element = typename LocalView::Element |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::FERequirementType = FERequirements_ |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::FlatBasis = typename Basis::FlatBasis |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::Geometry = typename Element::Geometry |
|
inherited |
|
inherited |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::GridView = typename FlatBasis::GridView |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis()) |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::LocalView = typename FlatBasis::LocalView |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::ResultRequirementsType = ResultRequirements<FERequirementType> |
|
inherited |
using Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::Traits = TraitsFromLocalView<LocalView, useEigenRef> |
|
inline |
VolumeLoad | The type for the volume load function. |
NeumannBoundaryLoad | The type for the Neumann boundary load function. |
globalBasis | The global basis for the finite element. |
element | The element for which the finite element is constructed. |
emod | Young's modulus. |
nu | Poisson's ratio. |
p_volumeLoad | Volume load function (default is LoadDefault). |
p_neumannBoundary | Neumann boundary patch (default is nullptr). |
p_neumannBoundaryLoad | Neumann boundary load function (default is LoadDefault). |
|
inline |
req | The ResultRequirementsType object specifying the requested results. |
local | Local position vector. |
result | Map to store the calculated results. |
|
inline |
par | The FERequirementType object. |
K | Matrix to store the calculated stiffness. |
|
inline |
par | The FERequirementType object. |
|
inlineprotected |
|
inline |
par | The FERequirementType object. |
force | Vector to store the calculated force. |
|
inlineprotected |
|
inline |
ScalarType | The scalar type for the displacement vector. |
par | The FERequirementType object. |
dx | Optional 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 |
|
inlineinherited |
|
inlineinherited |
|
inline |
|
inline |
par | The FERequirementType object. |
|
inlineconstexprinherited |
|
inline |
ScalarType | The scalar type for the strain vector. |
par | The FERequirementType object. |
dx | Optional displacement vector. |
double Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::emod_ |
std::shared_ptr<const Geometry> Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::geo_ |
Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType> > Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::localBasis |
|
staticconstexpr |
const BoundaryPatch<GridView>* Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::neumannBoundary |
std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)> Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::neumannBoundaryLoad |
double Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::nu_ |
|
staticconstexprinherited |
size_t Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::numberOfNodes {0} |
int Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::order {} |
std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)> Ikarus::LinearElastic< Basis_, FERequirements_, useEigenRef >::volumeLoad |
|
staticconstexpr |