Loading [MathJax]/extensions/tex2jax.js
OGS
MathLib::KelvinVector Namespace Reference

Detailed Description

The invariants and the Kelving mapping are explained in detail in the article "On Advantages of the Kelvin Mapping in Finite Element Implementations of Deformation Processes" [22].

Namespaces

 KelvinVector_detail
 

Classes

struct  Invariants
 

Typedefs

template<int DisplacementDim>
using KelvinVectorType = Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor >
 
template<int DisplacementDim>
using KelvinMatrixType = Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor >
 

Functions

template<>
Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > inverse (Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
 
template<>
Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > inverse (Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &v)
 
template<>
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor (Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
 
template<>
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor (Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &v)
 
template<>
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor (Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > const &v)
 
template<>
KelvinVectorType< 2 > tensorToKelvin< 2 > (Eigen::Matrix< double, 3, 3 > const &m)
 
template<>
KelvinVectorType< 3 > tensorToKelvin< 3 > (Eigen::Matrix< double, 3, 3 > const &m)
 
template<>
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor (Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
 
template<>
Eigen::Matrix< double, 6, 1 > kelvinVectorToSymmetricTensor (Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &v)
 
template<>
Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > kelvinVectorToSymmetricTensor (Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > const &v)
 
template<>
KelvinMatrixType< 2 > fourthOrderRotationMatrix< 2 > (Eigen::Matrix< double, 2, 2, Eigen::ColMajor, 2, 2 > const &transformation)
 
template<>
KelvinMatrixType< 3 > fourthOrderRotationMatrix< 3 > (Eigen::Matrix< double, 3, 3, Eigen::ColMajor, 3, 3 > const &transformation)
 
constexpr int kelvin_vector_dimensions (int const displacement_dim)
 Kelvin vector dimensions for given displacement dimension. More...
 
template<int KelvinVectorSize>
Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > inverse (Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &v)
 
template<int KelvinVectorSize>
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor (Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &v)
 
template<int DisplacementDim>
KelvinVectorType< DisplacementDim > tensorToKelvin (Eigen::Matrix< double, 3, 3 > const &m)
 
template<int KelvinVectorSize>
Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > kelvinVectorToSymmetricTensor (Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &v)
 
template<typename Derived >
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector (Eigen::MatrixBase< Derived > const &v)
 
template<int DisplacementDim>
KelvinVectorType< DisplacementDim > symmetricTensorToKelvinVector (std::vector< double > const &values)
 
template<int DisplacementDim>
KelvinMatrixType< DisplacementDim > fourthOrderRotationMatrix (Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::ColMajor, DisplacementDim, DisplacementDim > const &transformation)
 

Typedef Documentation

◆ KelvinMatrixType

template<int DisplacementDim>
using MathLib::KelvinVector::KelvinMatrixType = typedef Eigen::Matrix<double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor>

Kelvin matrix type for given displacement dimension.

Note
The Eigen matrix is always a fixed size vector in contrast to a shape matrix policy types like BMatrixPolicyType::KelvinMatrixType.

Definition at line 54 of file KelvinVector.h.

◆ KelvinVectorType

template<int DisplacementDim>
using MathLib::KelvinVector::KelvinVectorType = typedef Eigen::Matrix<double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor>

Kelvin vector type for given displacement dimension.

Note
The Eigen vector is always a fixed size vector in contrast to a shape matrix policy types like BMatrixPolicyType::KelvinVectorType.

Definition at line 46 of file KelvinVector.h.

Function Documentation

◆ fourthOrderRotationMatrix()

template<int DisplacementDim>
KelvinMatrixType<DisplacementDim> MathLib::KelvinVector::fourthOrderRotationMatrix ( Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::ColMajor, DisplacementDim, DisplacementDim > const &  transformation)

Rotation tensor for Kelvin mapped vectors and tensors. It is meant to be used for rotation of stress/strain tensors epsilon:Q and tangent stiffness tensors Q*C*Q^t. 2D and 3D implementations available.

Referenced by MaterialLib::Solids::LinearElasticOrthotropic< DisplacementDim >::getElasticTensor(), and MaterialLib::Solids::MFront::MFront< DisplacementDim >::integrateStress().

◆ fourthOrderRotationMatrix< 2 >()

template<>
KelvinMatrixType<2> MathLib::KelvinVector::fourthOrderRotationMatrix< 2 > ( Eigen::Matrix< double, 2, 2, Eigen::ColMajor, 2, 2 > const &  transformation)

Definition at line 162 of file KelvinVector.cpp.

186 {
187  // 1-based index access for convenience.
188  auto Q = [&](int const i, int const j)
189  { return transformation(i - 1, j - 1); };
190 
192  R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), 0,
193  std::sqrt(2) * Q(1, 1) * Q(1, 2), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
194  0, std::sqrt(2) * Q(2, 1) * Q(2, 2), 0, 0, 1, 0,
195  std::sqrt(2) * Q(1, 1) * Q(2, 1), std::sqrt(2) * Q(1, 2) * Q(2, 2), 0,
196  Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1);
197  return R;
198 }
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Definition: KelvinVector.h:56

References OGS_FATAL.

◆ fourthOrderRotationMatrix< 3 >()

template<>
KelvinMatrixType<3> MathLib::KelvinVector::fourthOrderRotationMatrix< 3 > ( Eigen::Matrix< double, 3, 3, Eigen::ColMajor, 3, 3 > const &  transformation)

Definition at line 162 of file KelvinVector.cpp.

203 {
204  // 1-based index access for convenience.
205  auto Q = [&](int const i, int const j)
206  { return transformation(i - 1, j - 1); };
207 
209  R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), Q(1, 3) * Q(1, 3),
210  std::sqrt(2) * Q(1, 1) * Q(1, 2), std::sqrt(2) * Q(1, 2) * Q(1, 3),
211  std::sqrt(2) * Q(1, 1) * Q(1, 3), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
212  Q(2, 3) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 2),
213  std::sqrt(2) * Q(2, 2) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 3),
214  Q(3, 1) * Q(3, 1), Q(3, 2) * Q(3, 2), Q(3, 3) * Q(3, 3),
215  std::sqrt(2) * Q(3, 1) * Q(3, 2), std::sqrt(2) * Q(3, 2) * Q(3, 3),
216  std::sqrt(2) * Q(3, 1) * Q(3, 3), std::sqrt(2) * Q(1, 1) * Q(2, 1),
217  std::sqrt(2) * Q(1, 2) * Q(2, 2), std::sqrt(2) * Q(1, 3) * Q(2, 3),
218  Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1),
219  Q(1, 2) * Q(2, 3) + Q(1, 3) * Q(2, 2),
220  Q(1, 1) * Q(2, 3) + Q(1, 3) * Q(2, 1), std::sqrt(2) * Q(2, 1) * Q(3, 1),
221  std::sqrt(2) * Q(2, 2) * Q(3, 2), std::sqrt(2) * Q(2, 3) * Q(3, 3),
222  Q(2, 1) * Q(3, 2) + Q(2, 2) * Q(3, 1),
223  Q(2, 2) * Q(3, 3) + Q(2, 3) * Q(3, 2),
224  Q(2, 1) * Q(3, 3) + Q(2, 3) * Q(3, 1), std::sqrt(2) * Q(1, 1) * Q(3, 1),
225  std::sqrt(2) * Q(1, 2) * Q(3, 2), std::sqrt(2) * Q(1, 3) * Q(3, 3),
226  Q(1, 1) * Q(3, 2) + Q(1, 2) * Q(3, 1),
227  Q(1, 2) * Q(3, 3) + Q(1, 3) * Q(3, 2),
228  Q(1, 1) * Q(3, 3) + Q(1, 3) * Q(3, 1);
229  return R;
230 }

◆ inverse() [1/3]

template<>
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> MathLib::KelvinVector::inverse ( Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &  v)

Definition at line 34 of file KelvinVector.cpp.

36 {
37  assert(Invariants<4>::determinant(v) != 0);
38 
39  Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> inv;
40  inv(0) = v(1) * v(2);
41  inv(1) = v(0) * v(2);
42  inv(2) = v(0) * v(1) - v(3) * v(3) / 2.;
43  inv(3) = -v(3) * v(2);
44  return inv / Invariants<4>::determinant(v);
45 }

References MathLib::KelvinVector::Invariants< KelvinVectorSize >::determinant().

Referenced by MaterialLib::Solids::Ehlers::calculatePlasticJacobian(), and MaterialLib::Solids::Ehlers::calculatePlasticResidual().

◆ inverse() [2/3]

template<>
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> MathLib::KelvinVector::inverse ( Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &  v)

Definition at line 48 of file KelvinVector.cpp.

50 {
51  assert(Invariants<6>::determinant(v) != 0);
52 
53  Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> inv;
54  inv(0) = v(1) * v(2) - v(4) * v(4) / 2.;
55  inv(1) = v(0) * v(2) - v(5) * v(5) / 2.;
56  inv(2) = v(0) * v(1) - v(3) * v(3) / 2.;
57  inv(3) = v(4) * v(5) / std::sqrt(2.) - v(3) * v(2);
58  inv(4) = v(3) * v(5) / std::sqrt(2.) - v(4) * v(0);
59  inv(5) = v(4) * v(3) / std::sqrt(2.) - v(1) * v(5);
60  return inv / Invariants<6>::determinant(v);
61 }

References MathLib::KelvinVector::Invariants< KelvinVectorSize >::determinant().

◆ inverse() [3/3]

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1> MathLib::KelvinVector::inverse ( Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &  v)

Inverse of a matrix in Kelvin vector representation. There are only implementations for the Kelvin vector size 4 and 6.

◆ kelvin_vector_dimensions()

constexpr int MathLib::KelvinVector::kelvin_vector_dimensions ( int const  displacement_dim)
constexpr

Kelvin vector dimensions for given displacement dimension.

Definition at line 23 of file KelvinVector.h.

24 {
25  if (displacement_dim == 2)
26  {
27  return 4;
28  }
29  else if (displacement_dim == 3)
30  {
31  return 6;
32  }
33  OGS_FATAL(
34  "Cannot convert displacement dimension {} to kelvin vector dimension.",
35  displacement_dim);
36 }
#define OGS_FATAL(...)
Definition: Error.h:26

References OGS_FATAL.

Referenced by ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::HydroMechanicsLocalAssembler(), ProcessLib::RichardsMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::IntegrationPointData(), ProcessLib::TH2M::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::IntegrationPointData(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::IntegrationPointData(), ProcessLib::ThermoRichardsMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::IntegrationPointData(), ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::PhaseFieldLocalAssembler(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::SmallDeformationLocalAssembler(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::SmallDeformationLocalAssemblerMatrix(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, IntegrationMethod, DisplacementDim >::SmallDeformationLocalAssemblerMatrixNearFracture(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::SmallDeformationNonlocalLocalAssembler(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ThermoMechanicalPhaseFieldLocalAssembler(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::ThermoMechanicsLocalAssembler(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assemble(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobianForDeformationEquations(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobianForHeatConductionEquations(), MaterialLib::Solids::Phasefield::calculateDegradedStress(), MaterialLib::Solids::Ehlers::calculateDResidualDEps(), MaterialLib::Solids::Ehlers::calculatePlasticJacobian(), MaterialLib::Solids::Ehlers::calculatePlasticResidual(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::RichardsMechanics::computeMicroPorosity(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::LIE::HydroMechanics::createHydroMechanicsProcess(), ProcessLib::HydroMechanics::createHydroMechanicsProcess(), ProcessLib::RichardsMechanics::createRichardsMechanicsProcess(), ProcessLib::SmallDeformation::createSmallDeformationProcess(), ProcessLib::TH2M::createTH2MProcess(), ProcessLib::ThermoHydroMechanics::createThermoHydroMechanicsProcess(), ProcessLib::ThermoMechanics::createThermoMechanicsProcess(), ProcessLib::ThermoRichardsMechanics::createThermoRichardsMechanicsProcess(), MaterialLib::Solids::MFront::MFront< DisplacementDim >::getBulkModulus(), MaterialLib::Solids::LinearElasticOrthotropic< DisplacementDim >::getElasticTensor(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getEpsilonMechanical(), ProcessLib::getIntegrationPointKelvinVectorData(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getIntPtSwellingStress(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSwellingStress(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getSigma(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getSigma(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getSigma(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::getSigma(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::getSigma(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getSwellingStress(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::getSwellingStress(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), MaterialLib::Solids::Ehlers::SolidEhlers< DisplacementDim >::integrateStress(), MaterialLib::Solids::MFront::MFront< DisplacementDim >::integrateStress(), MaterialLib::Solids::Ehlers::predict_sigma(), ProcessLib::setIntegrationPointKelvinVectorData(), symmetricTensorToKelvinVector(), MaterialLib::Solids::Lubby2::tangentStiffnessA(), and ProcessLib::RichardsMechanics::updateSwellingStressAndVolumetricStrain().

◆ kelvinVectorToSymmetricTensor() [1/4]

template<>
Eigen::Matrix<double, 4, 1> MathLib::KelvinVector::kelvinVectorToSymmetricTensor ( Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &  v)

Definition at line 142 of file KelvinVector.cpp.

144 {
145  Eigen::Matrix<double, 4, 1> m;
146  m << v[0], v[1], v[2], v[3] / std::sqrt(2.);
147  return m;
148 }

Referenced by ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assemble(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobian(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::assembleWithJacobianForPressureEquations(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getEpsilon(), ProcessLib::getIntegrationPointKelvinVectorData(), MaterialLib::Solids::Ehlers::SolidEhlers< DisplacementDim >::getInternalVariables(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::getIntPtSwellingStress(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSwellingStress(), and ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::updateConstitutiveVariables().

◆ kelvinVectorToSymmetricTensor() [2/4]

template<>
Eigen::Matrix<double, 6, 1> MathLib::KelvinVector::kelvinVectorToSymmetricTensor ( Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &  v)

Definition at line 151 of file KelvinVector.cpp.

153 {
154  Eigen::Matrix<double, 6, 1> m;
155  m << v[0], v[1], v[2], v[3] / std::sqrt(2.), v[4] / std::sqrt(2.),
156  v[5] / std::sqrt(2.);
157  return m;
158 }

◆ kelvinVectorToSymmetricTensor() [3/4]

template<>
Eigen::Matrix<double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1> MathLib::KelvinVector::kelvinVectorToSymmetricTensor ( Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > const &  v)

Definition at line 162 of file KelvinVector.cpp.

168 {
169  if (v.size() == 4)
170  {
171  return kelvinVectorToSymmetricTensor<4>(v);
172  }
173  if (v.size() == 6)
174  {
175  return kelvinVectorToSymmetricTensor<6>(v);
176  }
177  OGS_FATAL(
178  "Kelvin vector to tensor conversion expected an input vector of size 4 "
179  "or 6, but a vector of size {:d} was given.",
180  v.size());
181 }

◆ kelvinVectorToSymmetricTensor() [4/4]

template<int KelvinVectorSize>
Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1> MathLib::KelvinVector::kelvinVectorToSymmetricTensor ( Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &  v)

Conversion of a Kelvin vector to a short vector representation of a symmetric 3x3 matrix.

In the 2D case the entries for the xx, yy, zz, and xy components are stored. In the 3D case the entries for the xx, yy, zz, xy, yz, and xz components in that particular order are stored.

This is opposite of the symmetricTensorToKelvinVector()

Only implementations for KelvinVectorSize 4 and 6, and dynamic size vectors are provided.

◆ kelvinVectorToTensor() [1/4]

◆ kelvinVectorToTensor() [2/4]

template<>
Eigen::Matrix<double, 3, 3> MathLib::KelvinVector::kelvinVectorToTensor ( Eigen::Matrix< double, 6, 1, Eigen::ColMajor, 6, 1 > const &  v)

Definition at line 74 of file KelvinVector.cpp.

76 {
77  Eigen::Matrix<double, 3, 3> m;
78  m << v[0], v[3] / std::sqrt(2.), v[5] / std::sqrt(2.), v[3] / std::sqrt(2.),
79  v[1], v[4] / std::sqrt(2.), v[5] / std::sqrt(2.), v[4] / std::sqrt(2.),
80  v[2];
81  return m;
82 }

◆ kelvinVectorToTensor() [3/4]

template<>
Eigen::Matrix<double, 3, 3> MathLib::KelvinVector::kelvinVectorToTensor ( Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > const &  v)

Definition at line 85 of file KelvinVector.cpp.

91 {
92  if (v.size() == 4)
93  {
94  Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> v4;
95  v4 << v[0], v[1], v[2], v[3];
96  return kelvinVectorToTensor(v4);
97  }
98  if (v.size() == 6)
99  {
100  Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> v6;
101  v6 << v[0], v[1], v[2], v[3], v[4], v[5];
102  return kelvinVectorToTensor(v6);
103  }
104  OGS_FATAL(
105  "Conversion of dynamic Kelvin vector of size {:d} to a tensor is not "
106  "possible. Kelvin vector must be of size 4 or 6.",
107  v.size());
108 }
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor, Eigen::Dynamic, 1 > const &v)

References kelvinVectorToTensor(), and OGS_FATAL.

◆ kelvinVectorToTensor() [4/4]

template<int KelvinVectorSize>
Eigen::Matrix<double, 3, 3> MathLib::KelvinVector::kelvinVectorToTensor ( Eigen::Matrix< double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1 > const &  v)

Conversion of a Kelvin vector to a 3x3 matrix Only implementations for KelvinVectorSize 4 and 6 are provided.

◆ symmetricTensorToKelvinVector() [1/2]

template<typename Derived >
Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1> MathLib::KelvinVector::symmetricTensorToKelvinVector ( Eigen::MatrixBase< Derived > const &  v)

Conversion of a short vector representation of a symmetric 3x3 matrix to a Kelvin vector.

This is opposite of the kelvinVectorToSymmetricTensor()

Only implementations for KelvinVectorSize 4 and 6, and dynamic size vectors are provided.

Definition at line 170 of file KelvinVector.h.

171 {
172  static_assert(
173  (Eigen::MatrixBase<Derived>::ColsAtCompileTime == 1) ||
174  (Eigen::MatrixBase<Derived>::ColsAtCompileTime == Eigen::Dynamic),
175  "KelvinVector must be a column vector");
176  if (v.cols() != 1)
177  {
178  OGS_FATAL(
179  "KelvinVector must be a column vector, but input has {:d} columns.",
180  v.cols());
181  }
182 
183  Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1>
184  result;
185  if (v.rows() == 4)
186  {
187  result.resize(4, 1);
188  result << v[0], v[1], v[2], v[3] * std::sqrt(2.);
189  }
190  else if (v.rows() == 6)
191  {
192  result.resize(6, 1);
193  result << v[0], v[1], v[2], v[3] * std::sqrt(2.), v[4] * std::sqrt(2.),
194  v[5] * std::sqrt(2.);
195  }
196  else
197  {
198  OGS_FATAL(
199  "Symmetric tensor to Kelvin vector conversion expected an input "
200  "vector of size 4 or 6, but a vector of size {:d} was given.",
201  v.size());
202  }
203  return result;
204 }

References OGS_FATAL.

Referenced by ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, IntegrationMethod, DisplacementDim >::initializeConcrete(), ProcessLib::setIntegrationPointKelvinVectorData(), and symmetricTensorToKelvinVector().

◆ symmetricTensorToKelvinVector() [2/2]

template<int DisplacementDim>
KelvinVectorType<DisplacementDim> MathLib::KelvinVector::symmetricTensorToKelvinVector ( std::vector< double > const &  values)

Conversion of a short vector representation of a symmetric 3x3 matrix to a Kelvin vector.

This overload takes a std::vector for the tensor values.

Definition at line 211 of file KelvinVector.h.

213 {
214  constexpr int kelvin_vector_size =
215  kelvin_vector_dimensions(DisplacementDim);
216 
217  if (values.size() != kelvin_vector_size)
218  {
219  OGS_FATAL(
220  "Symmetric tensor to Kelvin vector conversion expected an input "
221  "vector of size {:d}, but a vector of size {:d} was given.",
222  kelvin_vector_size, values.size());
223  }
224 
226  Eigen::Map<typename MathLib::KelvinVector::KelvinVectorType<
227  DisplacementDim> const>(
228  values.data(), kelvin_vector_dimensions(DisplacementDim), 1));
229 }
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
KelvinVectorType< DisplacementDim > symmetricTensorToKelvinVector(std::vector< double > const &values)
Definition: KelvinVector.h:211
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23

References kelvin_vector_dimensions(), OGS_FATAL, and symmetricTensorToKelvinVector().

◆ tensorToKelvin()

template<int DisplacementDim>
KelvinVectorType<DisplacementDim> MathLib::KelvinVector::tensorToKelvin ( Eigen::Matrix< double, 3, 3 > const &  m)

Conversion of a 3x3 matrix to a Kelvin vector. Only implementations for KelvinVectorSize 4 and 6 are provided.

Referenced by ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, DisplacementDim >::updateConstitutiveVariables().

◆ tensorToKelvin< 2 >()

template<>
KelvinVectorType<2> MathLib::KelvinVector::tensorToKelvin< 2 > ( Eigen::Matrix< double, 3, 3 > const &  m)

Definition at line 85 of file KelvinVector.cpp.

112 {
113  assert(std::abs(m(0, 1) - m(1, 0)) <
114  std::numeric_limits<double>::epsilon());
115  assert(m(0, 2) == 0);
116  assert(m(1, 2) == 0);
117  assert(m(2, 0) == 0);
118  assert(m(2, 1) == 0);
119 
120  KelvinVectorType<2> v;
121  v << m(0, 0), m(1, 1), m(2, 2), m(0, 1) * std::sqrt(2.);
122  return v;
123 }

◆ tensorToKelvin< 3 >()

template<>
KelvinVectorType<3> MathLib::KelvinVector::tensorToKelvin< 3 > ( Eigen::Matrix< double, 3, 3 > const &  m)

Definition at line 85 of file KelvinVector.cpp.

127 {
128  assert(std::abs(m(0, 1) - m(1, 0)) <
129  std::numeric_limits<double>::epsilon());
130  assert(std::abs(m(1, 2) - m(2, 1)) <
131  std::numeric_limits<double>::epsilon());
132  assert(std::abs(m(0, 2) - m(2, 0)) <
133  std::numeric_limits<double>::epsilon());
134 
135  KelvinVectorType<3> v;
136  v << m(0, 0), m(1, 1), m(2, 2), m(0, 1) * std::sqrt(2.),
137  m(1, 2) * std::sqrt(2.), m(0, 2) * std::sqrt(2.);
138  return v;
139 }