KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | List of all members
HyperelasticMaterialLaw Class Reference

A hyperelastic material law. More...

#include <concepts.hh>

Detailed Description

A hyperelastic material law.

See also
Kaskade::MaterialLawBase

Definition at line 620 of file concepts.hh.

Public Types

using Scalar = unspecified
 The scalar field type (usually double). More...
 
using Tensor = unspecified
 The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim,dim>). More...
 
using VoigtTensor = unspecified
 The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim*(dim+1)/2,dim*(dim+1)/2>). More...
 

Public Member Functions

void setLinearizationPoint (Tensor const &e)
 Defines a new evaluation/linearization point. More...
 
Scalar d0 () const
 Evaluates the stored energy density \( W(\epsilon) \). More...
 
Scalar d1 (Tensor const &e1) const
 Evaluates the first directional derivative \( W'(\epsilon)\epsilon_1 \). More...
 
Scalar d2 (Tensor const &e1, Tensor const &e2) const
 Evaluates the second directional derivative \( W''(\epsilon)\epsilon_1\epsilon_2 \). More...
 
Tensor stress () const
 Returns the 2nd Piola-Kirchhoff stress tensor \( \sigma = W'(E) \) corresponding to the current strain \( E \). More...
 
VoigtTensor strainToStressMatrix () const
 Returns the tangent stiffness tensor \( C(E) = W''(E)\) mapping strain tensor variations to stress tensor (2nd Piola-Kirchhoff) variations. More...
 

Static Public Attributes

static int const dim = unspecified
 The spatial dimension (usually 2 or 3). More...
 

Member Typedef Documentation

◆ Scalar

using HyperelasticMaterialLaw::Scalar = unspecified

The scalar field type (usually double).

Definition at line 626 of file concepts.hh.

◆ Tensor

using HyperelasticMaterialLaw::Tensor = unspecified

The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim,dim>).

Definition at line 636 of file concepts.hh.

◆ VoigtTensor

The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim*(dim+1)/2,dim*(dim+1)/2>).

Definition at line 641 of file concepts.hh.

Member Function Documentation

◆ d0()

Scalar HyperelasticMaterialLaw::d0 ( ) const

Evaluates the stored energy density \( W(\epsilon) \).

◆ d1()

Scalar HyperelasticMaterialLaw::d1 ( Tensor const &  e1) const

Evaluates the first directional derivative \( W'(\epsilon)\epsilon_1 \).

◆ d2()

Scalar HyperelasticMaterialLaw::d2 ( Tensor const &  e1,
Tensor const &  e2 
) const

Evaluates the second directional derivative \( W''(\epsilon)\epsilon_1\epsilon_2 \).

The second derivative is symmetric, i.e. d2(e1,e2)==d2(e2,e1) holds for any e1, e2.

◆ setLinearizationPoint()

void HyperelasticMaterialLaw::setLinearizationPoint ( Tensor const &  e)

Defines a new evaluation/linearization point.

Parameters
ethe strain tensor

◆ strainToStressMatrix()

VoigtTensor HyperelasticMaterialLaw::strainToStressMatrix ( ) const

Returns the tangent stiffness tensor \( C(E) = W''(E)\) mapping strain tensor variations to stress tensor (2nd Piola-Kirchhoff) variations.

In this context, the symmetric matrices \( \sigma \) and \( \epsilon \) are interpreted as \(d(d+1)/2\)-vectors in Voigt notation as returned by Kaskade::pack. Note that strain tensor off-diagonal entries appear with factor two in the vector.

Now, a matrix \( C\in\mathbb{R}^{(d+1)d/2\times(d+1)d/2} \) is returned such that for infinitesimal changes \( \delta\epsilon \) the corresponding change in \( \sigma(\epsilon+\delta\epsilon) = \sigma(\epsilon)+\delta\sigma\) is given as \( \delta\sigma = C \delta\epsilon \). Note that due to the symmetry of the stored energy Hessian \( W''(E) \), the Voigt representation of the stiffness tensor is a symmetric matrix.

The following invariant holds for any strain tensors eps1 and eps2:

d2(eps1,eps2) == pack(eps1)*(strainToStressMatrix()*pack(eps2))
VoigtTensor strainToStressMatrix() const
Returns the tangent stiffness tensor mapping strain tensor variations to stress tensor (2nd Piola-Ki...
Scalar d2(Tensor const &e1, Tensor const &e2) const
Evaluates the second directional derivative .
Dune::FieldVector< Scalar, n *(n+1)/2 > pack(Dune::FieldMatrix< Scalar, n, n > const &e, Scalar s=2.0)
Returns the Voigt vector representation of a symmetric matrix (with or without doubled off-diagonal e...

The Drucker/Hill stability criterion states that the returned matrix shall be positive semidefinite.

◆ stress()

Tensor HyperelasticMaterialLaw::stress ( ) const

Returns the 2nd Piola-Kirchhoff stress tensor \( \sigma = W'(E) \) corresponding to the current strain \( E \).

The following invariant holds for any strain tensor eps:

d1(eps) == contraction(stress(),eps)
Scalar d1(Tensor const &e1) const
Evaluates the first directional derivative .
Tensor stress() const
Returns the 2nd Piola-Kirchhoff stress tensor corresponding to the current strain .
T contraction(Dune::FieldMatrix< T, n, m > const &A, Dune::FieldMatrix< T, n, m > const &B)
Matrix contraction (Frobenius product) .
Definition: fixdune.hh:494

Member Data Documentation

◆ dim

int const HyperelasticMaterialLaw::dim = unspecified
static

The spatial dimension (usually 2 or 3).

Definition at line 631 of file concepts.hh.


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