KASKADE 7 development version
Public Member Functions | List of all members
Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim > Class Template Reference

Convenience class for handling orthotropic materials in linear elastomechanics. More...

#include <elastoVariationalFunctionals.hh>

Detailed Description

template<class Scalar, int dim>
class Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >

Convenience class for handling orthotropic materials in linear elastomechanics.

This class defines the variational formulation of the linear elasticity for orthotropic materials for vector-valued displacements \( u \), i.e. the variational functional \( u \mapsto \frac{1}{2} \epsilon : \mathcal{C}\epsilon \) with \( \epsilon = \frac{1}{2}(u_x + u_x^T) \). In an appropriately rotated coordinate system, the stiffness matrix, representing the tensor \(\mathcal{C}\) by interpreting \( \epsilon, \sigma \) as vectors, has the form

\[ C = \begin{bmatrix} * & * & \\ * & * & \\ & & * \end{bmatrix}, \quad C = \begin{bmatrix} * & * & * & & & \\ * & * & * & & & \\ * & * & * & & & \\ & & & * & & \\ & & & & * & \\ & & & & & * \end{bmatrix} \]

depending on the spatial dimension (see Wikipedia).

Definition at line 429 of file elastoVariationalFunctionals.hh.

Public Member Functions

 OrthotropicLameNavier (Dune::FieldMatrix< Scalar, dim, dim > const &orth_, Dune::FieldMatrix< Scalar, dim, dim > const &mat_)
 Constructor. More...
 
 OrthotropicLameNavier (Dune::FieldMatrix< Scalar, dim, dim > const &orth_, Dune::FieldMatrix< Scalar, dim, dim > const &mat1, Dune::FieldVector< Scalar, dim *(dim-1)/2 > const &mat2)
 Constructor. More...
 
 OrthotropicLameNavier (ElasticModulus const &p)
 Default constructor. More...
 
 OrthotropicLameNavier ()
 Default constructor. More...
 
Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *(dim+1)/2 > const & stiffnessMatrix () const
 Returns the stiffness tensor. More...
 
void setLinearizationPoint (Dune::FieldMatrix< Scalar, dim, dim > const &du0_)
 Defines the displacement derivative around which to linearize. More...
 
Scalar d0 () const
 Computes the elastic energy \( \frac{1}{2} \epsilon(u_x) : \mathcal{C}\epsilon(u_x) \). More...
 
Dune::FieldVector< Scalar, dim > d1 (VariationalArg< Scalar, dim > const &arg) const
 Computes the elastic energy \( \frac{1}{2} \epsilon(u_x) : \mathcal{C}\epsilon(u_x) \). More...
 
Dune::FieldMatrix< Scalar, dim, dim > d2 (VariationalArg< Scalar, dim > const &arg1, VariationalArg< Scalar, dim > const &arg2) const
 Computes the local Hessian of the variational functional. 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. More...
 
void setLinear (bool lin)
 

Constructor & Destructor Documentation

◆ OrthotropicLameNavier() [1/4]

template<class Scalar , int dim>
Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::OrthotropicLameNavier ( Dune::FieldMatrix< Scalar, dim, dim > const &  orth_,
Dune::FieldMatrix< Scalar, dim, dim > const &  mat_ 
)

Constructor.

Parameters
orthThe material coordinate system matrix. Has to be orthonormal.
matThe material parameters.

The columns of the orth matrix are the material coordinate directions, i.e. a multiplication by this matrix transforms from material coordinates to world coordinates.

mat is a matrix containing the material parameters in the following way:

\[ \begin{bmatrix} E_1 & G_{12} \\ \nu_{12} & E_2 \end{bmatrix}, \quad \begin{bmatrix} E_1 & G_{12} & G_{13} \\ \nu_{12} & E_2 & G_{23} \\ \nu_{13} & \nu_{23} & E_3 \end{bmatrix}, \]

where \( E_i \) is the elastic modulus in material coordinate direction \( i \), \( \nu_{ij} \) is the Poissons ratio for contraction in \( j \) direction due to strain in \( i \) direction, and \( G_{ij} \) is the shear modulus in \( ij \) plane.

◆ OrthotropicLameNavier() [2/4]

template<class Scalar , int dim>
Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::OrthotropicLameNavier ( Dune::FieldMatrix< Scalar, dim, dim > const &  orth_,
Dune::FieldMatrix< Scalar, dim, dim > const &  mat1,
Dune::FieldVector< Scalar, dim *(dim-1)/2 > const &  mat2 
)

Constructor.

Parameters
orthThe material coordinate system matrix. Has to be orthonormal.
mat1material parameters
mat2material parameters

The columns of the orth matrix are the material coordinate directions, i.e. a multiplication by this matrix transforms from material coordinates to world coordinates.

mat1 is a \( d \times d \) matrix containing the top left part of the stiffness matrix \( C \), and mat2 is a \( d(d-1)/2 \) vector containing the diagonal of the bottom right part.

\[ C = \begin{bmatrix} \mathrm{mat1} & \\ & \mathrm{diag}(\mathrm{mat2}) \end{bmatrix}, \]

◆ OrthotropicLameNavier() [3/4]

template<class Scalar , int dim>
Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::OrthotropicLameNavier ( ElasticModulus const &  p)

Default constructor.

This initializes the material to an isotropic material with given material parameters.

◆ OrthotropicLameNavier() [4/4]

template<class Scalar , int dim>
Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::OrthotropicLameNavier ( )
inline

Default constructor.

This initializes the material to an isotropic material with Lame parameters set to 1.

Definition at line 476 of file elastoVariationalFunctionals.hh.

Member Function Documentation

◆ d0()

template<class Scalar , int dim>
Scalar Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::d0 ( ) const

Computes the elastic energy \( \frac{1}{2} \epsilon(u_x) : \mathcal{C}\epsilon(u_x) \).

◆ d1()

template<class Scalar , int dim>
Dune::FieldVector< Scalar, dim > Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::d1 ( VariationalArg< Scalar, dim > const &  arg) const

Computes the elastic energy \( \frac{1}{2} \epsilon(u_x) : \mathcal{C}\epsilon(u_x) \).

arg is assumed to be a scalar variational argument.

◆ d2()

template<class Scalar , int dim>
Dune::FieldMatrix< Scalar, dim, dim > Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::d2 ( VariationalArg< Scalar, dim > const &  arg1,
VariationalArg< Scalar, dim > const &  arg2 
) const

Computes the local Hessian of the variational functional. 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.

◆ setLinear()

template<class Scalar , int dim>
void Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::setLinear ( bool  lin)
inline

Definition at line 520 of file elastoVariationalFunctionals.hh.

◆ setLinearizationPoint()

template<class Scalar , int dim>
void Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::setLinearizationPoint ( Dune::FieldMatrix< Scalar, dim, dim > const &  du0_)
inline

Defines the displacement derivative around which to linearize.

Definition at line 492 of file elastoVariationalFunctionals.hh.

◆ stiffnessMatrix()

template<class Scalar , int dim>
Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *(dim+1)/2 > const & Kaskade::Elastomechanics::OrthotropicLameNavier< Scalar, dim >::stiffnessMatrix ( ) const
inline

Returns the stiffness tensor.

The stiffness tensor \( C \) is represented in matrix form for Voigt notation (i.e. the mapping from Voigt strain vector to Voigt stress vector). In 3D this reads

\[ \begin{bmatrix} \sigma_{11} \\ \sigma_{22} \\ \sigma_{33} \\ \sigma_{23} \\ \sigma_{13} \\ \sigma_{12} \end{bmatrix} = C \begin{bmatrix} \epsilon_{11} \\ \epsilon_{22} \\ \epsilon_{33} \\ 2\epsilon_{23} \\ 2\epsilon_{13} \\ 2\epsilon_{12} \end{bmatrix} \]

Definition at line 487 of file elastoVariationalFunctionals.hh.


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