Namespaces | Enumerations | Functions | Variables
Marmot::ContinuumMechanics::VoigtNotation Namespace Reference

Namespaces

 Derivatives
 
 Invariants
 
 Transformations
 

Enumerations

enum  VoigtSize { OneD = 1, TwoD = 3, ThreeD = 6, Axial = 4 }
 

Functions

constexpr VoigtSize voigtSizeFromDimension (int x)
 
Eigen::Vector3d voigtToPlaneVoigt (const Marmot::Vector6d &voigt)
 
Marmot::Vector6d planeVoigtToVoigt (const Eigen::Vector3d &voigtPlane)
 
template<enum VoigtSize voigtSize>
Eigen::Matrix< double, voigtSize, 1 > reduce3DVoigt (const Marmot::Vector6d &Voigt3D)
 
template<enum VoigtSize voigtSize>
Marmot::Vector6d make3DVoigt (const Eigen::Matrix< double, voigtSize, 1 > &Voigt)
 
template<typename T >
Eigen::Matrix< T, 3, 3 > voigtToStrain (const Eigen::Matrix< T, 6, 1 > &voigt)
 
template<typename T >
Eigen::Matrix< T, 3, 3 > voigtToStress (const Eigen::Matrix< T, 6, 1 > &voigt)
 
Marmot::Vector6d strainToVoigt (const Eigen::Matrix3d &strainTensor)
 
template<typename T = double>
Eigen::Matrix< T, 6, 1 > stressToVoigt (const Eigen::Matrix< T, 3, 3 > &stressTensor)
 
template<int nDim>
Eigen::Matrix< double, nDim, nDim > stressMatrixFromVoigt (const Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > &Voigt)
 
template<int nDim>
Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > voigtFromStrainMatrix (const Eigen::Matrix< double, nDim, nDim > &strain)
 
Vector6d strainToVoigt (const Matrix3d &strainTensor)
 
Vector6d planeVoigtToVoigt (const Vector3d &voigtPlane)
 

Variables

const Marmot::Vector6d P = ( Vector6d() << 1, 1, 1, 2, 2, 2 ).finished()
 
const Marmot::Vector6d PInv = ( Vector6d() << 1, 1, 1, .5, .5, .5 ).finished()
 
const Marmot::Vector6d I = ( Vector6d() << 1, 1, 1, 0, 0, 0 ).finished()
 
const Marmot::Vector6d IHyd = ( Vector6d() << 1. / 3, 1. / 3, 1. / 3, 0, 0, 0 ).finished()
 
const Matrix6d IDev
 

Enumeration Type Documentation

◆ VoigtSize

Enumerator
OneD 
TwoD 
ThreeD 
Axial 

Function Documentation

◆ voigtSizeFromDimension()

constexpr VoigtSize Marmot::ContinuumMechanics::VoigtNotation::voigtSizeFromDimension ( int  x)
constexpr

◆ voigtToPlaneVoigt()

Vector3d Marmot::ContinuumMechanics::VoigtNotation::voigtToPlaneVoigt ( const Marmot::Vector6d voigt)

Converts a 3D voigt notated vector to a plane stress vector.

\[ \displaystyle \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ \sigma_{33} = 0\\ \sigma_{12}\\ \sigma_{13} = 0\\ \sigma_{23} = 0 \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ \sigma_{12} \end{bmatrix} \]

◆ planeVoigtToVoigt() [1/2]

Marmot::Vector6d Marmot::ContinuumMechanics::VoigtNotation::planeVoigtToVoigt ( const Eigen::Vector3d &  voigtPlane)

Converts a voigt notated plane stress vector to a 3D vector.

\[ \displaystyle \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ \sigma_{12} \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ 0\\ \sigma_{12}\\ 0\\ 0 \end{bmatrix} \]

◆ reduce3DVoigt()

template<enum VoigtSize voigtSize>
Eigen::Matrix< double, voigtSize, 1 > Marmot::ContinuumMechanics::VoigtNotation::reduce3DVoigt ( const Marmot::Vector6d Voigt3D)

Reduces a 3D voigt notated vector to a lower dimension defined by the template parameter 'voigtSize'.

Considered cases:

  • voigtSize \( = 1 \):

    \[ \displaystyle \begin{bmatrix} \sigma_{11}\\ \sigma_{22} = 0\\ \sigma_{33} = 0\\ \sigma_{12} = 0\\ \sigma_{13} = 0\\ \sigma_{23} = 0 \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \sigma_{11} \]

  • voigtSize \( = 3 \): Calls function voigtToPlaneVoigt().
  • voigtSize \( = 6 \): Returns the input vector

◆ make3DVoigt()

template<enum VoigtSize voigtSize>
Marmot::Vector6d Marmot::ContinuumMechanics::VoigtNotation::make3DVoigt ( const Eigen::Matrix< double, voigtSize, 1 > &  Voigt)

Computes a 3D voigt notated vector from a vector of lower dimension defined by the template parameter 'voigtSize'.

Considered cases:

  • voigtSize \( = 1 \):

    \[ \displaystyle \sigma_{11} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \sigma_{11}\\ 0\\ 0\\ 0\\ 0\\ 0 \end{bmatrix} \]

  • voigtSize \( = 3 \): Calls function planeVoigtToVoigt().
  • voigtSize \( = 6 \): Returns the input vector

◆ voigtToStrain()

template<typename T >
Eigen::Matrix< T, 3, 3 > Marmot::ContinuumMechanics::VoigtNotation::voigtToStrain ( const Eigen::Matrix< T, 6, 1 > &  voigt)

Converts a voigt notated strain vector to its corresponding strain tensor

\[ \displaystyle \begin{bmatrix} \varepsilon_{11}\\ \varepsilon_{22}\\ \varepsilon_{33}\\ \gamma_{12}\\ \gamma_{13}\\ \gamma_{23} \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \varepsilon_{11} & \gamma_{12}/2 & \gamma_{13}/2 \\ \gamma_{12}/2 & \varepsilon_{22} & \gamma_{23}/2 \\ \gamma_{13}/2 & \gamma_{23}/2 & \varepsilon_{33} \\ \end{bmatrix} \]

◆ voigtToStress()

template<typename T >
Eigen::Matrix< T, 3, 3 > Marmot::ContinuumMechanics::VoigtNotation::voigtToStress ( const Eigen::Matrix< T, 6, 1 > &  voigt)

Converts a voigt notated stress vector to its corresponding stress tensor

\[ \displaystyle \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ \sigma_{33}\\ \sigma_{12}\\ \sigma_{13}\\ \sigma_{23} \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \sigma_{11} & \sigma_{12} & \sigma_{13} \\ \sigma_{12} & \sigma_{22} & \sigma_{23} \\ \sigma_{13} & \sigma_{23} & \sigma_{33} \\ \end{bmatrix} \]

◆ strainToVoigt() [1/2]

Marmot::Vector6d Marmot::ContinuumMechanics::VoigtNotation::strainToVoigt ( const Eigen::Matrix3d &  strainTensor)

Converts a strain tensor to its corresponding voigt notated strain vector

\[ \displaystyle \begin{bmatrix} \varepsilon_{11} & \varepsilon_{12} & \varepsilon_{13}\\ \varepsilon_{21} & \varepsilon_{22} & \varepsilon_{23} \\ \varepsilon_{31} & \varepsilon_{32} & \varepsilon_{33} \\ \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \varepsilon_{11}\\ \varepsilon_{22}\\ \varepsilon_{33}\\ 2\,\varepsilon_{12}\\ 2\,\varepsilon_{13}\\ 2\,\varepsilon_{23} \end{bmatrix} \]

◆ stressToVoigt()

template<typename T = double>
Eigen::Matrix< T, 6, 1 > Marmot::ContinuumMechanics::VoigtNotation::stressToVoigt ( const Eigen::Matrix< T, 3, 3 > &  stressTensor)

Converts a stress tensor to its corresponding voigt notated stress vector

\[ \displaystyle \begin{bmatrix} \sigma_{11} & \sigma_{12} & \sigma_{13}\\ \sigma_{21} & \sigma_{22} & \sigma_{23} \\ \sigma_{31} & \sigma_{32} & \sigma_{33} \\ \end{bmatrix} \hspace{.5cm} \Rightarrow \hspace{.5cm} \begin{bmatrix} \sigma_{11}\\ \sigma_{22}\\ \sigma_{33}\\ \sigma_{12}\\ \sigma_{13}\\ \sigma_{23} \end{bmatrix} \]

◆ stressMatrixFromVoigt()

template<int nDim>
Eigen::Matrix< double, nDim, nDim > Marmot::ContinuumMechanics::VoigtNotation::stressMatrixFromVoigt ( const Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > &  Voigt)

◆ voigtFromStrainMatrix()

template<int nDim>
Eigen::Matrix< double, VOIGTFROMDIM( nDim ), 1 > Marmot::ContinuumMechanics::VoigtNotation::voigtFromStrainMatrix ( const Eigen::Matrix< double, nDim, nDim > &  strain)

◆ strainToVoigt() [2/2]

Vector6d Marmot::ContinuumMechanics::VoigtNotation::strainToVoigt ( const Matrix3d strainTensor)

◆ planeVoigtToVoigt() [2/2]

Vector6d Marmot::ContinuumMechanics::VoigtNotation::planeVoigtToVoigt ( const Vector3d voigtPlane)

Variable Documentation

◆ P

const Vector6d Marmot::ContinuumMechanics::VoigtNotation::P = ( Vector6d() << 1, 1, 1, 2, 2, 2 ).finished()

◆ PInv

const Vector6d Marmot::ContinuumMechanics::VoigtNotation::PInv = ( Vector6d() << 1, 1, 1, .5, .5, .5 ).finished()

◆ I

const Vector6d Marmot::ContinuumMechanics::VoigtNotation::I = ( Vector6d() << 1, 1, 1, 0, 0, 0 ).finished()

◆ IHyd

const Vector6d Marmot::ContinuumMechanics::VoigtNotation::IHyd = ( Vector6d() << 1. / 3, 1. / 3, 1. / 3, 0, 0, 0 ).finished()

◆ IDev

const Matrix6d Marmot::ContinuumMechanics::VoigtNotation::IDev
Initial value:
= ( Matrix6d() <<
2./3, -1./3, -1./3, 0, 0, 0,
-1./3, 2./3, -1./3, 0, 0, 0,
-1./3, -1./3, 2./3, 0, 0, 0,
0, 0, 0, 1, 0, 0,
0, 0, 0, 0, 1, 0,
0, 0, 0, 0, 0, 1).finished()
Marmot::Matrix6d
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: MarmotTypedefs.h:35