OGS
MathLib::KelvinVector::Invariants< KelvinVectorSize > Struct Template Referencefinal

Detailed Description

template<int KelvinVectorSize>
struct MathLib::KelvinVector::Invariants< KelvinVectorSize >

Invariants used in mechanics, based on Kelvin representation of the vectors and matrices. The invariants are computed at process creation time.

Definition at line 87 of file KelvinVector.h.

#include <KelvinVector.h>

Public Member Functions

double determinant (Eigen::Matrix< double, 6, 1 > const &v)
double determinant (Eigen::Matrix< double, 4, 1 > const &v)

Static Public Member Functions

static double determinant (Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
 Determinant of a matrix in Kelvin vector representation.
static double equivalentStress (Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
static double FrobeniusNorm (Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
 Get the norm of the deviatoric stress.
static double J2 (Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
static double J3 (Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
static double trace (Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
 Trace of the corresponding tensor.
static Eigen::Vector3d diagonal (Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)

Static Public Attributes

static Eigen::Matrix< double, KelvinVectorSize, KelvinVectorSize > const deviatoric_projection
static Eigen::Matrix< double, KelvinVectorSize, KelvinVectorSize > const spherical_projection
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
 Kelvin mapping of 2nd order identity tensor.
static Eigen::Matrix< double, KelvinVectorSize, 1 > const ones2

Member Function Documentation

◆ determinant() [1/3]

double MathLib::KelvinVector::Invariants< 4 >::determinant ( Eigen::Matrix< double, 4, 1 > const & v)

Definition at line 21 of file KelvinVector.cpp.

22{
23 return v(2) * (v(0) * v(1) - v(3) * v(3) / 2.);
24}
static const double v

References MathLib::v.

◆ determinant() [2/3]

double MathLib::KelvinVector::Invariants< 6 >::determinant ( Eigen::Matrix< double, 6, 1 > const & v)

Definition at line 13 of file KelvinVector.cpp.

14{
15 return v(0) * v(1) * v(2) + v(3) * v(4) * v(5) / std::sqrt(2.) -
16 v(3) * v(3) * v(2) / 2. - v(4) * v(4) * v(0) / 2. -
17 v(5) * v(5) * v(1) / 2.;
18}

References MathLib::v.

◆ determinant() [3/3]

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::determinant ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & v)
static

Determinant of a matrix in Kelvin vector representation.

References MathLib::v.

Referenced by MathLib::KelvinVector::inverse(), MathLib::KelvinVector::inverse(), and J3().

◆ diagonal()

template<int KelvinVectorSize>
Eigen::Vector3d MathLib::KelvinVector::Invariants< KelvinVectorSize >::diagonal ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & v)
static

Diagonal of the corresponding tensor which is always of length 3 in 2D and 3D cases.

Definition at line 47 of file KelvinVector-impl.h.

49{
50 return v.template topLeftCorner<3, 1>();
51}

References MathLib::v.

Referenced by equivalentStress(), J2(), J3(), and trace().

◆ equivalentStress()

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::equivalentStress ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & deviatoric_v)
static

The von Mises equivalent stress.

Note
The input vector must have trace equal zero.

Definition at line 11 of file KelvinVector-impl.h.

13{
15 4e-12 * diagonal(deviatoric_v).norm());
16 return std::sqrt(3 * J2(deviatoric_v));
17}
static Eigen::Vector3d diagonal(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
static double trace(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Trace of the corresponding tensor.
static double J2(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)

References diagonal(), J2(), and trace().

◆ FrobeniusNorm()

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::FrobeniusNorm ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & deviatoric_v)
static

◆ J2()

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::J2 ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & deviatoric_v)
static

Second invariant of deviatoric tensor.

Note
The input vector must have trace equal zero.

Definition at line 27 of file KelvinVector-impl.h.

29{
31 4e-12 * diagonal(deviatoric_v).norm());
32 return 0.5 * deviatoric_v.transpose() * deviatoric_v;
33}

References diagonal(), and trace().

Referenced by equivalentStress().

◆ J3()

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::J3 ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & deviatoric_v)
static

Third invariant, equal to determinant of a deviatoric tensor.

Note
The input vector must have trace equal zero.

Definition at line 38 of file KelvinVector-impl.h.

40{
42 4e-12 * diagonal(deviatoric_v).norm());
44}
static double determinant(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Determinant of a matrix in Kelvin vector representation.

References determinant(), diagonal(), and trace().

◆ trace()

template<int KelvinVectorSize>
double MathLib::KelvinVector::Invariants< KelvinVectorSize >::trace ( Eigen::Matrix< double, KelvinVectorSize, 1 > const & v)
static

Trace of the corresponding tensor.

Definition at line 55 of file KelvinVector-impl.h.

57{
58 return diagonal(v).sum();
59}

References diagonal(), and MathLib::v.

Referenced by ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianEvalConstitutiveSetting(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianForDeformationEquations(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForPressureEquations(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianHydroEquations(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianPhaseFieldEquations(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), equivalentStress(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity(), J2(), J3(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcrete(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), and ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::updateConstitutiveRelations().

Member Data Documentation

◆ deviatoric_projection

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, KelvinVectorSize> const MathLib::KelvinVector::Invariants< KelvinVectorSize >::deviatoric_projection
static

Kelvin mapping of deviatoric projection tensor. \(A_{\rm dev} = P_{\rm dev}:A\) for \(A\) being a second order tensor.

Definition at line 95 of file KelvinVector.h.

Referenced by ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianHydroEquations().

◆ identity2

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, 1> const MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2
static

Kelvin mapping of 2nd order identity tensor.

Definition at line 101 of file KelvinVector.h.

Referenced by ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleBlockMatricesWithJacobian(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianForDeformationEquations(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForPressureEquations(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianHydroEquations(), ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianPhaseFieldEquations(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::convertInitialStressType(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), and ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::updateConstitutiveRelations().

◆ ones2

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, 1> const MathLib::KelvinVector::Invariants< KelvinVectorSize >::ones2
static

Kelvin mapping of a 2nd order tensor whose elements are one. In 3D case, it is equal to \([1., 1., 1., \sqrt(2), \sqrt(2), \sqrt(2)]\).

Definition at line 105 of file KelvinVector.h.

Referenced by ProcessLib::HMPhaseField::HMPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianHydroEquations().

◆ spherical_projection

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, KelvinVectorSize> const MathLib::KelvinVector::Invariants< KelvinVectorSize >::spherical_projection
static

The documentation for this struct was generated from the following files: