KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | Protected Attributes | List of all members
Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor > Class Template Reference

General base class for variational functionals defined in terms of hyperelastic stored energies. More...

#include <elastoVariationalFunctionals.hh>

Detailed Description

template<class HyperelasticEnergy, class StrainTensor>
class Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >

General base class for variational functionals defined in terms of hyperelastic stored energies.

This defines the interface for stored energies working on displacement derivatives \( W = W(u_x) \). The interface is intended to be used for defining the domain cache of variational functionals.

The class has as a member an object of the hyperelastic energy defining the stored energy \( W(E) \) depending on the Green-Lagrange strain tensor.

Template Parameters
HyperelasticEnergythe hyperelastic energy class
StrainTensorthe strain tensor class, usually one of GreenLagrangeTensor or LinearizedGreenLagrangeTensor

Example for linear elastomechanics of steel:

energy.setLinearizationPoint(du);
static ElasticModulus const & material(std::string const &name)
Returns the material parameters of the given materials.
General base class for variational functionals defined in terms of hyperelastic stored energies.
Linearized (right) Green-Lagrange strain tensor, the workhorse of linear elastomechanics.
Definition: elasto.hh:280

Definition at line 43 of file elastoVariationalFunctionals.hh.

Public Types

using Scalar = typename HyperelasticEnergy::Scalar
 
using Tensor = Dune::FieldMatrix< Scalar, dim, dim >
 

Public Member Functions

template<typename... Args>
 HyperelasticVariationalFunctional (Args &&... args)
 Constructor. More...
 
 HyperelasticVariationalFunctional (HyperelasticEnergy const &energy_, StrainTensor const &strain_)
 Constructor. More...
 
void setLinearizationPoint (Tensor const &du0)
 Defines the displacement derivative around which to linearize. More...
 
Scalar d0 () const
 Computes the hyperelastic stored energy density. More...
 
Dune::FieldVector< Scalar, dimd1 (VariationalArg< Scalar, dim > const &arg) const
 Computes the directional derivative of the hyperelastic stored energy density. More...
 
Tensor d2 (VariationalArg< Scalar, dim > const &arg1, VariationalArg< Scalar, dim > const &arg2) const
 Computes the second directional derivative of the hyperelastic stored energy density. More...
 
HyperelasticEnergy & getHyperelasticEnergy ()
 Provides access to the stored energy density. More...
 
HyperelasticEnergy const & getHyperelasticEnergy () const
 getHyperelasticEnergy provides read access to the stored energy density. More...
 
Tensor const & getDu () const
 getDu provides read access to the current displacement derivative. More...
 
StrainTensor const & getStrain () const
 getStrain provides read access to the current strain tensor. More...
 
Tensor stress () const
 Returns the 2nd Piola-Kirchhoff stress tensor \( S = \frac{\partial W}{\partial E} \) corresponding to the current displacement derivative. More...
 
Tensor cauchyStress () const
 Returns the Cauchy stress tensor \( \sigma = \frac{1}{J} F S F^T \) corresponding to the current displacement derivative. More...
 

Static Public Attributes

static int const dim = HyperelasticEnergy::dim
 

Protected Attributes

Tensor du
 
HyperelasticEnergy energy
 
StrainTensor strain
 

Member Typedef Documentation

◆ Scalar

template<class HyperelasticEnergy , class StrainTensor >
using Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::Scalar = typename HyperelasticEnergy::Scalar

Definition at line 46 of file elastoVariationalFunctionals.hh.

◆ Tensor

template<class HyperelasticEnergy , class StrainTensor >
using Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::Tensor = Dune::FieldMatrix<Scalar,dim,dim>

Definition at line 48 of file elastoVariationalFunctionals.hh.

Constructor & Destructor Documentation

◆ HyperelasticVariationalFunctional() [1/2]

template<class HyperelasticEnergy , class StrainTensor >
template<typename... Args>
Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::HyperelasticVariationalFunctional ( Args &&...  args)
inline

Constructor.

The constructor arguments are forwarded to the hyperelastic stored energy constructor. The strain tensor is default initialized.

Definition at line 57 of file elastoVariationalFunctionals.hh.

◆ HyperelasticVariationalFunctional() [2/2]

template<class HyperelasticEnergy , class StrainTensor >
Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::HyperelasticVariationalFunctional ( HyperelasticEnergy const &  energy_,
StrainTensor const &  strain_ 
)
inline

Constructor.

Parameters
energya hyperelastic energy functional, will be copied to an internal energy member variable
straina strain tensor functional, will be copied to an internal strain variable

Definition at line 67 of file elastoVariationalFunctionals.hh.

Member Function Documentation

◆ cauchyStress()

template<class HyperelasticEnergy , class StrainTensor >
Tensor Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::cauchyStress ( ) const
inline

Returns the Cauchy stress tensor \( \sigma = \frac{1}{J} F S F^T \) corresponding to the current displacement derivative.

Here, \( S \) is the second Piola-Kirchhoff stress tensor and \( F = I + u_x \) the deformation derivative.

The Cauchy stress \( \sigma \), also known as "true" stress, describes the tractions in the deformed coordinates, and therefore the "true" force density on a (cutting) surface.

Definition at line 207 of file elastoVariationalFunctionals.hh.

◆ d0()

template<class HyperelasticEnergy , class StrainTensor >
Scalar Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::d0 ( ) const
inline

Computes the hyperelastic stored energy density.

Definition at line 91 of file elastoVariationalFunctionals.hh.

◆ d1()

template<class HyperelasticEnergy , class StrainTensor >
Dune::FieldVector< Scalar, dim > Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::d1 ( VariationalArg< Scalar, dim > const &  arg) const
inline

Computes the directional derivative of the hyperelastic stored energy density.

Definition at line 99 of file elastoVariationalFunctionals.hh.

◆ d2()

template<class HyperelasticEnergy , class StrainTensor >
Tensor Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::d2 ( VariationalArg< Scalar, dim > const &  arg1,
VariationalArg< Scalar, dim > const &  arg2 
) const
inline

Computes the second directional derivative of the hyperelastic stored energy density.

For scalar functions \( \phi, \psi \), the values and derivatives at a certain point are provided in arg1 and arg 2, this computes the \( d\times d \) matrix \( A \) with \( A_{ij} = \epsilon((u_i)_x)^T \mathcal{C} \epsilon((v_j)_x) \), where \( u_i = \phi e_i, v_j = \psi e_j \) and \( e_k \) is the \( k \)-th unit vector.

This is exactly what is to be returned from the d2 method of a a variational functional's domain cache implementing the Lame-Navier equations.

Definition at line 127 of file elastoVariationalFunctionals.hh.

◆ getDu()

template<class HyperelasticEnergy , class StrainTensor >
Tensor const & Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::getDu ( ) const
inline

getDu provides read access to the current displacement derivative.

Definition at line 177 of file elastoVariationalFunctionals.hh.

Referenced by Kaskade::Elastomechanics::FirstPiolaKirchhoffStress< HyperelasticEnergy, StrainTensor >::deformationDerivative().

◆ getHyperelasticEnergy() [1/2]

template<class HyperelasticEnergy , class StrainTensor >
HyperelasticEnergy & Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::getHyperelasticEnergy ( )
inline

Provides access to the stored energy density.

The stored energy function has the linearization (and evaluation) point given by the latest call to the setLinearizationPoint() method.

Definition at line 161 of file elastoVariationalFunctionals.hh.

Referenced by Kaskade::Elastomechanics::FirstPiolaKirchhoffStress< HyperelasticEnergy, StrainTensor >::energyD2De().

◆ getHyperelasticEnergy() [2/2]

template<class HyperelasticEnergy , class StrainTensor >
HyperelasticEnergy const & Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::getHyperelasticEnergy ( ) const
inline

getHyperelasticEnergy provides read access to the stored energy density.

Definition at line 169 of file elastoVariationalFunctionals.hh.

◆ getStrain()

template<class HyperelasticEnergy , class StrainTensor >
StrainTensor const & Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::getStrain ( ) const
inline

◆ setLinearizationPoint()

template<class HyperelasticEnergy , class StrainTensor >
void Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::setLinearizationPoint ( Tensor const &  du0)
inline

Defines the displacement derivative around which to linearize.

This method shall be called when the evaluation of a stored energy and its derivatives is intended for a different displacement derivative, in particular in the "evaluateAt" method of a domain cache.

Parameters
du0the current displacement derivative \( u_x \)

Definition at line 81 of file elastoVariationalFunctionals.hh.

◆ stress()

template<class HyperelasticEnergy , class StrainTensor >
Tensor Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::stress ( ) const
inline

Returns the 2nd Piola-Kirchhoff stress tensor \( S = \frac{\partial W}{\partial E} \) corresponding to the current displacement derivative.

This is the stress in the reference (material) coordinates, and does not provide the actual force density in the deformed configuration. See cauchyStress.

Definition at line 197 of file elastoVariationalFunctionals.hh.

Referenced by Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::cauchyStress(), Kaskade::Elastomechanics::FirstPiolaKirchhoffStress< HyperelasticEnergy, StrainTensor >::d0(), and Kaskade::Elastomechanics::FirstPiolaKirchhoffStress< HyperelasticEnergy, StrainTensor >::d1().

Member Data Documentation

◆ dim

template<class HyperelasticEnergy , class StrainTensor >
int const Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::dim = HyperelasticEnergy::dim
static

◆ du

template<class HyperelasticEnergy , class StrainTensor >
Tensor Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::du
protected

◆ energy

template<class HyperelasticEnergy , class StrainTensor >
HyperelasticEnergy Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::energy
protected

◆ strain

template<class HyperelasticEnergy , class StrainTensor >
StrainTensor Kaskade::Elastomechanics::HyperelasticVariationalFunctional< HyperelasticEnergy, StrainTensor >::strain
protected

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