KASKADE 7 development version
|
Convenience class for handling linear elastomechanics. More...
#include <elastoVariationalFunctionals.hh>
Convenience class for handling linear elastomechanics.
This class defines the variational formulation of the Lame-Navier operator \( u \mapsto -2\mu\Delta u - \lambda \nabla \mbox{div} u\) for vector-valued displacements \( u \), i.e. the variational functional \( u \mapsto \frac{\lambda}{2} \mbox{tr}(E)^2 + \mu(E:E) \) with \( E = \frac{1}{2}(u_x + u_x^T) \).
Definition at line 234 of file elastoVariationalFunctionals.hh.
Public Types | |
using | Scalar = Scalar_ |
using | Tensor = Dune::FieldMatrix< Scalar, dim, dim > |
Public Member Functions | |
LameNavier () | |
Constructor. Both Lame constants are initialized to 1, the linearization point to \( (u_0)_x = 0 \). More... | |
LameNavier (ElasticModulus const &moduli) | |
Dune::FieldVector< Scalar, dim > | d1 (VariationalArg< Scalar, dim > const &arg) const |
d1 As the base method, computes the directional derivative of the hyperelastic stored energy density. More... | |
Base::Tensor | d2 (VariationalArg< Scalar, dim > const &arg1, VariationalArg< Scalar, dim > const &arg2) const |
d2 As the base method, computes the second directional derivative of the hyperelastic stored energy density. But this transformed version (only valid in the completely linear case) is more efficient. More... | |
void | setElasticModuli (ElasticModulus const &modulus) |
Define material properties. More... | |
ElasticModulus | getElasticModuli () |
Retrieve current material properties. More... | |
Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *(dim+1)/2 > | strainToStressMatrix () const |
Returns the matrix \( C \) mapping the strain tensor to the stress tensor: \( \sigma = C \epsilon \). More... | |
Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *dim > | displacementDerivativeToStrainMatrix () const |
returns the matrix \( E \) mapping the displacement derivative to the strain tensor: \( \epsilon = E u_x \) 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, dim > | d1 (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... | |
MaterialLaws::StVenantKirchhoff< dim_, double > & | getHyperelasticEnergy () |
Provides access to the stored energy density. More... | |
MaterialLaws::StVenantKirchhoff< dim_, double > 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... | |
LinearizedGreenLagrangeTensor< double, dim_ > 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 Member Functions | |
static Dune::FieldVector< Scalar, dim > | mv (Dune::FieldVector< Scalar, dim *(dim+1)/2 > const &A, Dune::FieldVector< Scalar, dim > const &x) |
Computes the symmetric matrix times vector product, when the matrix is compactly stored in a dim*(dim+1)/2-vector. More... | |
Static Public Attributes | |
static int const | dim = dim_ |
Protected Attributes | |
Tensor | du |
MaterialLaws::StVenantKirchhoff< dim_, double > | energy |
LinearizedGreenLagrangeTensor< double, dim_ > | strain |
using Kaskade::Elastomechanics::LameNavier< dim_, Scalar_ >::Scalar = Scalar_ |
Definition at line 242 of file elastoVariationalFunctionals.hh.
|
inherited |
Definition at line 48 of file elastoVariationalFunctionals.hh.
|
inline |
Constructor. Both Lame constants are initialized to 1, the linearization point to \( (u_0)_x = 0 \).
Definition at line 254 of file elastoVariationalFunctionals.hh.
|
inline |
Definition at line 257 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
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.
|
inlineinherited |
Computes the hyperelastic stored energy density.
Definition at line 91 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
Computes the directional derivative of the hyperelastic stored energy density.
Definition at line 99 of file elastoVariationalFunctionals.hh.
|
inline |
d1 As the base method, computes the directional derivative of the hyperelastic stored energy density.
But this transformed version (only valid in the completely linear case) is probably more efficient.
Definition at line 267 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
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.
|
inline |
d2 As the base method, computes the second directional derivative of the hyperelastic stored energy density. But this transformed version (only valid in the completely linear case) is more efficient.
Definition at line 280 of file elastoVariationalFunctionals.hh.
|
inline |
returns the matrix \( E \) mapping the displacement derivative to the strain tensor: \( \epsilon = E u_x \)
In this context, the symmetric matrix \( \epsilon \) is interpreted as vector in \( \mathbb{R}^{d(d+1)/2} \):
\[ \epsilon = \begin{bmatrix} \epsilon_{11} \\ \epsilon_{22} \\ \epsilon_{33} \\ \epsilon_{12} \\ \epsilon_{13} \\ \epsilon_{23} \end{bmatrix} \quad\mbox{or}\quad \epsilon = \begin{bmatrix} \epsilon_{11} \\ \epsilon_{22} \\ \epsilon_{12} \end{bmatrix} \]
The displacement derivative \( u_x \) is interpreted as vector in \( \mathbb{R}^{d^2} \) (column-major format):
\[ \frac{\partial u_i}{\partial x_j} = (u_x)_{i+dj} \]
Definition at line 387 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
getDu provides read access to the current displacement derivative.
Definition at line 177 of file elastoVariationalFunctionals.hh.
|
inline |
Retrieve current material properties.
Definition at line 316 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
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.
|
inlineinherited |
getHyperelasticEnergy provides read access to the stored energy density.
Definition at line 169 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
getStrain provides read access to the current strain tensor.
Definition at line 185 of file elastoVariationalFunctionals.hh.
|
inlinestatic |
Computes the symmetric matrix times vector product, when the matrix is compactly stored in a dim*(dim+1)/2-vector.
This is useful, e.g., for computing normal stresses \( \sigma n \).
Definition at line 355 of file elastoVariationalFunctionals.hh.
|
inline |
Define material properties.
This invalidates previously set linearization points.
Definition at line 306 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
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.
du0 | the current displacement derivative \( u_x \) |
Definition at line 81 of file elastoVariationalFunctionals.hh.
|
inline |
Returns the matrix \( C \) mapping the strain tensor to the stress tensor: \( \sigma = C \epsilon \).
In this contex, the symmetric matrices \( \sigma \) and \( \epsilon \) are interpreted as vectors in \( \mathbb{R}^{d(d+1)/2} \):
\[ \sigma = \begin{bmatrix} \sigma_{11} \\ \sigma_{22} \\ \sigma_{33} \\ \sigma_{12} \\ \sigma_{13} \\ \sigma_{23} \end{bmatrix} \quad\mbox{or}\quad \sigma = \begin{bmatrix} \sigma_{11} \\ \sigma_{22} \\ \sigma_{12} \end{bmatrix} \]
Definition at line 328 of file elastoVariationalFunctionals.hh.
|
inlineinherited |
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.
|
static |
Definition at line 248 of file elastoVariationalFunctionals.hh.
Referenced by Kaskade::Elastomechanics::LameNavier< dim_, Scalar_ >::d2(), Kaskade::Elastomechanics::LameNavier< dim_, Scalar_ >::displacementDerivativeToStrainMatrix(), Kaskade::Elastomechanics::LameNavier< dim_, Scalar_ >::mv(), and Kaskade::Elastomechanics::LameNavier< dim_, Scalar_ >::strainToStressMatrix().
|
protectedinherited |
Definition at line 217 of file elastoVariationalFunctionals.hh.
|
protectedinherited |
Definition at line 218 of file elastoVariationalFunctionals.hh.
|
protectedinherited |
Definition at line 219 of file elastoVariationalFunctionals.hh.