version 0.4.1
Ikarus::NonLinearElastic< PreFE, FE, PRE > Class Template Reference

NonLinearElastic class represents a non-linear elastic finite element. More...

#include <ikarus/finiteelements/mechanics/nonlinearelastic.hh>

Public Types

using Traits = PreFE::Traits
 
using Basis = typename Traits::Basis
 
using FlatBasis = typename Traits::FlatBasis
 
using FERequirementType = typename Traits::FERequirementType
 
using LocalView = typename Traits::LocalView
 
using Geometry = typename Traits::Geometry
 
using GridView = typename Traits::GridView
 
using Element = typename Traits::Element
 
using Material = PRE::Material
 
using Pre = PRE
 
using LocalBasisType = decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis())
 
using SupportedResultTypes = std::tuple< decltype(makeRT< ResultTypes::PK2Stress >())>
 

Public Member Functions

 NonLinearElastic (const Pre &pre)
 Constructor for the NonLinearElastic class. More...
 
template<typename ScalarType = double>
auto displacementFunction (const FERequirementType &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
 Get the displacement function for the given FERequirementType. More...
 
template<typename ScalarType = double>
auto strainFunction (const FERequirementType &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
 The strain function for the given FERequirementType. More...
 
template<typename ScalarType , int strainDim, bool voigt = true>
auto materialTangent (const Eigen::Vector< ScalarType, strainDim > &strain) const
 Get the material tangent for the given strain. More...
 
template<typename ScalarType , int strainDim>
auto getInternalEnergy (const Eigen::Vector< ScalarType, strainDim > &strain) const
 Get the internal energy for the given strain. More...
 
template<typename ScalarType , int strainDim, bool voigt = true>
auto getStress (const Eigen::Vector< ScalarType, strainDim > &strain) const
 Get the stress for the given strain. More...
 
const Geometrygeometry () const
 
size_t numberOfNodes () const
 
int order () const
 
template<template< typename, int, int > class RT>
requires (canProvideResultType<RT>())
auto calculateAtImpl (const FERequirementType &req, const Dune::FieldVector< double, Traits::mydim > &local, Dune::PriorityTag< 1 >) const
 Calculates a requested result at a specific local position. More...
 

Static Public Attributes

static constexpr int myDim = Traits::mydim
 
static constexpr auto strainType = StrainTags::greenLagrangian
 

Protected Member Functions

void bindImpl ()
 A helper function to bind the local view to the element. More...
 
template<typename ScalarType >
void calculateMatrixImpl (const FERequirementType &par, typename Traits::template MatrixType<> K, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
 Calculate the matrix associated with the given FERequirementType. More...
 
template<typename ScalarType >
auto calculateScalarImpl (const FERequirementType &par, const std::optional< std::reference_wrapper< 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< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
 

Detailed Description

template<typename PreFE, typename FE, typename PRE>
class Ikarus::NonLinearElastic< PreFE, FE, PRE >
Template Parameters
PreFEThe type of the total pre finite element.
FEThe type of the finite element.
PREThe type of the non-linear elastic pre finite element.

Member Typedef Documentation

◆ Basis

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Basis = typename Traits::Basis

◆ Element

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Element = typename Traits::Element

◆ FERequirementType

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::FERequirementType = typename Traits::FERequirementType

◆ FlatBasis

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::FlatBasis = typename Traits::FlatBasis

◆ Geometry

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Geometry = typename Traits::Geometry

◆ GridView

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::GridView = typename Traits::GridView

◆ LocalBasisType

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis())

◆ LocalView

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::LocalView = typename Traits::LocalView

◆ Material

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Material = PRE::Material

◆ Pre

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Pre = PRE

◆ SupportedResultTypes

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::SupportedResultTypes = std::tuple<decltype(makeRT<ResultTypes::PK2Stress>())>

◆ Traits

template<typename PreFE , typename FE , typename PRE >
using Ikarus::NonLinearElastic< PreFE, FE, PRE >::Traits = PreFE::Traits

Constructor & Destructor Documentation

◆ NonLinearElastic()

template<typename PreFE , typename FE , typename PRE >
Ikarus::NonLinearElastic< PreFE, FE, PRE >::NonLinearElastic ( const Pre pre)
inlineexplicit
Parameters
preThe pre fe

Member Function Documentation

◆ bindImpl()

template<typename PreFE , typename FE , typename PRE >
void Ikarus::NonLinearElastic< PreFE, FE, PRE >::bindImpl ( )
inlineprotected

◆ calculateAtImpl()

template<typename PreFE , typename FE , typename PRE >
template<template< typename, int, int > class RT>
requires (canProvideResultType<RT>())
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::calculateAtImpl ( const FERequirementType req,
const Dune::FieldVector< double, Traits::mydim > &  local,
Dune::PriorityTag< 1 >   
) const
inline
Parameters
reqThe FERequirementType object holding the global solution.
localLocal position vector.
Returns
calculated result
Template Parameters
RTThe type representing the requested result.

◆ calculateMatrixImpl()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType >
void Ikarus::NonLinearElastic< PreFE, FE, PRE >::calculateMatrixImpl ( const FERequirementType par,
typename Traits::template MatrixType<>  K,
const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &  dx = std::nullopt 
) const
inlineprotected
Template Parameters
ScalarTypeThe scalar type for the calculation.
Parameters
parThe FERequirementType object specifying the requirements for the calculation.
KThe matrix to store the calculated result.

◆ calculateScalarImpl()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType >
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::calculateScalarImpl ( const FERequirementType par,
const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &  dx = std::nullopt 
) const -> ScalarType
inlineprotected

◆ calculateVectorImpl()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType >
void Ikarus::NonLinearElastic< PreFE, FE, PRE >::calculateVectorImpl ( const FERequirementType par,
typename Traits::template VectorType< ScalarType >  force,
const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &  dx = std::nullopt 
) const
inlineprotected

◆ displacementFunction()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType = double>
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::displacementFunction ( const FERequirementType par,
const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &  dx = std::nullopt 
) const
inline
Template Parameters
ScalarTypeThe scalar type for the displacement function.
Parameters
parThe FERequirementType object.
dxOptional displacement vector.
Returns
A StandardLocalFunction representing the displacement function.
Here is the caller graph for this function:

◆ geometry()

template<typename PreFE , typename FE , typename PRE >
const Geometry & Ikarus::NonLinearElastic< PreFE, FE, PRE >::geometry ( ) const
inline

◆ getInternalEnergy()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType , int strainDim>
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::getInternalEnergy ( const Eigen::Vector< ScalarType, strainDim > &  strain) const
inline
Template Parameters
ScalarTypeThe scalar type for the material and strain.
strainDimThe dimension of the strain vector.
Parameters
strainThe strain vector.
Returns
The internal energy calculated using the material's storedEnergy function.
Here is the caller graph for this function:

◆ getStress()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType , int strainDim, bool voigt = true>
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::getStress ( const Eigen::Vector< ScalarType, strainDim > &  strain) const
inline
Template Parameters
ScalarTypeThe scalar type for the material and strain.
strainDimThe dimension of the strain vector.
voigtA boolean indicating whether to use the Voigt notation for stress.
Parameters
strainThe strain vector.
Returns
The stress vector calculated using the material's stresses function.
Here is the caller graph for this function:

◆ materialTangent()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType , int strainDim, bool voigt = true>
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::materialTangent ( const Eigen::Vector< ScalarType, strainDim > &  strain) const
inline
Template Parameters
ScalarTypeThe scalar type for the material and strain.
strainDimThe dimension of the strain vector.
voigtFlag indicating whether to use Voigt notation.
Parameters
strainThe strain vector.
Returns
The material tangent calculated using the material's tangentModuli function.
Here is the caller graph for this function:

◆ numberOfNodes()

template<typename PreFE , typename FE , typename PRE >
size_t Ikarus::NonLinearElastic< PreFE, FE, PRE >::numberOfNodes ( ) const
inline

◆ order()

template<typename PreFE , typename FE , typename PRE >
int Ikarus::NonLinearElastic< PreFE, FE, PRE >::order ( ) const
inline

◆ strainFunction()

template<typename PreFE , typename FE , typename PRE >
template<typename ScalarType = double>
auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::strainFunction ( const FERequirementType par,
const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &  dx = std::nullopt 
) const
inline
Template Parameters
ScalarTypeThe scalar type for the strain function.
Parameters
parThe FERequirementType object.
dxOptional displacement vector.
Returns
The strain function calculated using greenLagrangeStrains.
Here is the caller graph for this function:

Member Data Documentation

◆ myDim

template<typename PreFE , typename FE , typename PRE >
constexpr int Ikarus::NonLinearElastic< PreFE, FE, PRE >::myDim = Traits::mydim
staticconstexpr

◆ strainType

template<typename PreFE , typename FE , typename PRE >
constexpr auto Ikarus::NonLinearElastic< PreFE, FE, PRE >::strainType = StrainTags::greenLagrangian
staticconstexpr

The documentation for this class was generated from the following file: