KASKADE 7 development version
Modules | Namespaces | Classes | Functions

Predefined standard differential operators and material laws to be used as building blocks for elastomechanics problems. More...

Modules

 Stationary elasticity
 Building blocks for stationary elasticity problems (linear elasticity and hyperelasticity)
 
 Viscoplasticity
 Building blocks for viscoplastic problems.
 
 Contact Problems
 Components for defining elastomechanical contact problems.
 

Namespaces

namespace  Kaskade::Elastomechanics::MaterialLaws
 A namespace containing various material laws.
 

Classes

class  Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >
 A lower-dimensional kinematics ansatz useful for shell models. More...
 
class  Kaskade::Elastomechanics::ShellKinematics< dim, Real >
 A lower-dimensional kinematics ansatz useful for shell models. More...
 
class  Kaskade::Elastomechanics::ShellEnergy< Energy, Real >
 A stored energy density class for Taylor-based shell models. More...
 

Functions

template<class Scalar , int n>
Dune::FieldVector< Scalar, n *(n+1)/2 > Kaskade::Elastomechanics::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 entries). More...
 
template<class Scalar , int n>
Dune::FieldMatrix< Scalar, ElastoDetails::dim(n), ElastoDetails::dim(n)> Kaskade::Elastomechanics::unpack (Dune::FieldVector< Scalar, n > const &c, Scalar s=2.0)
 Returns the symmetric matrix representation of a Voigt vector with or without doubled off-diagonal entries. More...
 
template<int d, class Scalar >
Dune::FieldMatrix< Scalar, d, d > Kaskade::Elastomechanics::deviatoricPart (Dune::FieldMatrix< Scalar, d, d > const &s)
 The deviatoric part of a tensor. More...
 

Detailed Description

Predefined standard differential operators and material laws to be used as building blocks for elastomechanics problems.

Function Documentation

◆ deviatoricPart()

template<int d, class Scalar >
Dune::FieldMatrix< Scalar, d, d > Kaskade::Elastomechanics::deviatoricPart ( Dune::FieldMatrix< Scalar, d, d > const &  s)

The deviatoric part of a tensor.

This computes the deviatoric part \( \mathrm{dev}(s) = s - \frac{I_1}{d} I \) of the tensor \( s \), where \( I_1 = \mathrm{tr}(s) \) is the first invariant.

Definition at line 173 of file elasto.hh.

◆ pack()

template<class Scalar , int n>
Dune::FieldVector< Scalar, n *(n+1)/2 > Kaskade::Elastomechanics::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 entries).

The result is

\[ \begin{bmatrix} a_{11} & a_{12} \\ * & a_{22} \end{bmatrix} \rightarrow \begin{bmatrix} a_{11} \\ a_{22} \\ s a_{12} \end{bmatrix}, \qquad \begin{bmatrix} a_{11} & a_{12} & a_{13} \\ * & a_{22} & a_{23} \\ * & * & a_{33} \end{bmatrix} \rightarrow \begin{bmatrix} a_{11} \\ a_{22} \\ a_{33} \\ sa_{23} \\ s a_{13} \\ sa_{12} \end{bmatrix} \]

Note that traditional Voigt notation treats the doubling of off-diagonal entries differently for stress and strain tensors ( \( s=2\) for strains but \( s=1 \) for stresses, such that \(\sigma:\epsilon = \sigma^T \epsilon\)).

The inverse operation is provided as unpack.

Parameters
ethe symmetric matrix
sthe off-diagonal factor (defaults to 2)

Referenced by Kaskade::Elastomechanics::MaterialLaws::OrthotropicLinearMaterial< dim, Scalar >::d0(), Kaskade::Elastomechanics::MaterialLaws::OrthotropicLinearMaterial< dim, Scalar >::d1(), and Kaskade::Elastomechanics::MaterialLaws::OrthotropicLinearMaterial< dim, Scalar >::d2().

◆ unpack()

template<class Scalar , int n>
Dune::FieldMatrix< Scalar, ElastoDetails::dim(n), ElastoDetails::dim(n)> Kaskade::Elastomechanics::unpack ( Dune::FieldVector< Scalar, n > const &  c,
Scalar  s = 2.0 
)

Returns the symmetric matrix representation of a Voigt vector with or without doubled off-diagonal entries.

This is the inverse operation of pack.

Parameters
ethe symmetric matrix
sthe off-diagonal factor (defaults to 2)

Referenced by Kaskade::Elastomechanics::MaterialLaws::MaterialLawBase< dim_, Scalar_, MaterialLaw >::strainToStressMatrix().