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" [30].

Namespaces

namespace  KelvinVector_detail
 

Classes

struct  Invariants
 

Typedefs

template<int DisplacementDim>
using KelvinVectorType
 
template<int DisplacementDim>
using KelvinMatrixType
 

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<>
Eigen::Matrix< double, 2, 4 > liftVectorToKelvin< 2 > (Eigen::Matrix< double, 2, 1 > const &v)
 
template<>
Eigen::Matrix< double, 3, 6 > liftVectorToKelvin< 3 > (Eigen::Matrix< double, 3, 1 > const &v)
 
template<>
Eigen::Matrix< double, 2, 1 > reduceKelvinToVector< 2 > (Eigen::Matrix< double, 2, 4 > const &m)
 
template<>
Eigen::Matrix< double, 3, 1 > reduceKelvinToVector< 3 > (Eigen::Matrix< double, 3, 6 > const &m)
 
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.
 
template<int DisplacementDim>
constexpr auto KVzero ()
 Returns an expressions for a Kelvin vector filled with zero.
 
template<int DisplacementDim>
constexpr auto KMzero ()
 Returns an expressions for a Kelvin matrix filled with zero.
 
template<int DisplacementDim>
constexpr auto KVnan ()
 Returns an expressions for a Kelvin vector filled with NaN.
 
template<int DisplacementDim>
constexpr auto KMnan ()
 Returns an expressions for a Kelvin matrix filled with NaN.
 
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>
Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> liftVectorToKelvin (Eigen::Matrix< double, DisplacementDim, 1 > const &v)
 
template<int DisplacementDim>
Eigen::Matrix< double, DisplacementDim, 1 > reduceKelvinToVector (Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> const &m)
 
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
Initial value:
Eigen::Matrix<double, kelvin_vector_dimensions(DisplacementDim),
kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor>
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.

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 55 of file KelvinVector.h.

◆ KelvinVectorType

template<int DisplacementDim>
using MathLib::KelvinVector::KelvinVectorType
Initial value:
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 47 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.

◆ 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.

243{
244 // 1-based index access for convenience.
245 auto Q = [&](int const i, int const j)
246 { return transformation(i - 1, j - 1); };
247
249 R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), 0,
250 std::sqrt(2) * Q(1, 1) * Q(1, 2), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
251 0, std::sqrt(2) * Q(2, 1) * Q(2, 2), 0, 0, 1, 0,
252 std::sqrt(2) * Q(1, 1) * Q(2, 1), std::sqrt(2) * Q(1, 2) * Q(2, 2), 0,
253 Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1);
254 return R;
255}
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType

References kelvinVectorToSymmetricTensor(), OGS_FATAL, and MathLib::v.

◆ 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.

260{
261 // 1-based index access for convenience.
262 auto Q = [&](int const i, int const j)
263 { return transformation(i - 1, j - 1); };
264
266 R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), Q(1, 3) * Q(1, 3),
267 std::sqrt(2) * Q(1, 1) * Q(1, 2), std::sqrt(2) * Q(1, 2) * Q(1, 3),
268 std::sqrt(2) * Q(1, 1) * Q(1, 3), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
269 Q(2, 3) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 2),
270 std::sqrt(2) * Q(2, 2) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 3),
271 Q(3, 1) * Q(3, 1), Q(3, 2) * Q(3, 2), Q(3, 3) * Q(3, 3),
272 std::sqrt(2) * Q(3, 1) * Q(3, 2), std::sqrt(2) * Q(3, 2) * Q(3, 3),
273 std::sqrt(2) * Q(3, 1) * Q(3, 3), std::sqrt(2) * Q(1, 1) * Q(2, 1),
274 std::sqrt(2) * Q(1, 2) * Q(2, 2), std::sqrt(2) * Q(1, 3) * Q(2, 3),
275 Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1),
276 Q(1, 2) * Q(2, 3) + Q(1, 3) * Q(2, 2),
277 Q(1, 1) * Q(2, 3) + Q(1, 3) * Q(2, 1), std::sqrt(2) * Q(2, 1) * Q(3, 1),
278 std::sqrt(2) * Q(2, 2) * Q(3, 2), std::sqrt(2) * Q(2, 3) * Q(3, 3),
279 Q(2, 1) * Q(3, 2) + Q(2, 2) * Q(3, 1),
280 Q(2, 2) * Q(3, 3) + Q(2, 3) * Q(3, 2),
281 Q(2, 1) * Q(3, 3) + Q(2, 3) * Q(3, 1), std::sqrt(2) * Q(1, 1) * Q(3, 1),
282 std::sqrt(2) * Q(1, 2) * Q(3, 2), std::sqrt(2) * Q(1, 3) * Q(3, 3),
283 Q(1, 1) * Q(3, 2) + Q(1, 2) * Q(3, 1),
284 Q(1, 2) * Q(3, 3) + Q(1, 3) * Q(3, 2),
285 Q(1, 1) * Q(3, 3) + Q(1, 3) * Q(3, 1);
286 return R;
287}

◆ 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(), and MathLib::v.

Referenced by MaterialLib::Solids::Ehlers::thetaSigmaDerivatives().

◆ 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(), and MathLib::v.

◆ 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()

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

Kelvin vector dimensions for given displacement dimension.

Definition at line 24 of file KelvinVector.h.

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

References OGS_FATAL.

Referenced by ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssembler(), ProcessLib::ThermoHydroMechanics::IntegrationPointData< BMatricesType, ShapeMatrixTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, NPoints >::IntegrationPointData(), ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::PhaseFieldLocalAssembler(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerMatrix(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerMatrixNearFracture(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicalPhaseFieldLocalAssembler(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::ThermoMechanicsLocalAssembler(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianEvalConstitutiveSetting(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForDeformationEquations(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobianForHeatConductionEquations(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::assembleWithJacobianSingleIP(), MaterialLib::Solids::Ehlers::calculateDResidualDEps(), MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStress(), MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStressWithRankineEnergy(), MaterialLib::Solids::Phasefield::calculateOrthoVolDevDegradedStress(), MaterialLib::Solids::Ehlers::calculatePlasticJacobian(), MaterialLib::Solids::Ehlers::calculatePlasticResidual(), MaterialLib::Solids::Phasefield::calculateVolDevDegradedStress(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), ProcessLib::RichardsMechanics::computeMicroPorosity(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::LIE::HydroMechanics::createHydroMechanicsProcess(), ProcessLib::createInitialStress(), ProcessLib::LargeDeformation::createLargeDeformationProcess(), ProcessLib::SmallDeformation::createSmallDeformationProcess(), MaterialPropertyLib::Function::Implementation< D >::createSymbolTable(), ProcessLib::ThermoMechanics::createThermoMechanicsProcess(), ProcessLib::TH2M::ConstitutiveRelations::PorosityModelNonConstantSolidPhaseVolumeFraction< DisplacementDim >::dEval(), ProcessLib::TH2M::ConstitutiveRelations::SolidDensityModelNonConstantSolidPhaseVolumeFraction< DisplacementDim >::dEval(), ProcessLib::TH2M::ConstitutiveRelations::PermeabilityModel< DisplacementDim >::eval(), ProcessLib::TH2M::ConstitutiveRelations::PorosityModelNonConstantSolidPhaseVolumeFraction< DisplacementDim >::eval(), ProcessLib::TH2M::ConstitutiveRelations::SolidDensityModelNonConstantSolidPhaseVolumeFraction< DisplacementDim >::eval(), ProcessLib::TH2M::ConstitutiveRelations::SolidThermalExpansionModel< DisplacementDim >::eval(), ProcessLib::TH2M::ConstitutiveRelations::TotalStressModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature::SolidMechanicsModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature::SwellingModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::FluidThermalExpansionModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::PermeabilityModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::PorosityModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::TransportPorosityModel< DisplacementDim >::eval(), MaterialLib::Solids::MFront::MFrontGeneric< DisplacementDim, Gradients, TDynForces, ExtStateVars >::getBulkModulus(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilon(), ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getEpsilon(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilon(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getEpsilon(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilonM(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getEpsilonMechanical(), ProcessLib::getIntegrationPointKelvinVectorData(), ProcessLib::getIntegrationPointKelvinVectorData(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getSigma(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::getSigma(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::Reflection::detail::GetFlattenedIPDataFromLocAsm< Dim, Accessor_IPDataVecInLocAsm, Accessor_CurrentLevelFromIPDataVecElement >::operator()(), MaterialLib::Solids::Ehlers::predict_sigma(), MaterialLib::Solids::MFront::Variable< Derived >::rows(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), ProcessLib::setIntegrationPointKelvinVectorData(), ProcessLib::Reflection::detail::setIPData(), MaterialLib::Solids::MFront::OGSMFrontTangentOperatorBlocksView< DisplacementDim, ForcesGradsCombinations >::size(), symmetricTensorToKelvinVector(), MaterialLib::Solids::MFront::tangentOperatorDataMFrontToOGS(), MaterialLib::Solids::Lubby2::tangentStiffnessA(), MaterialLib::Solids::Ehlers::thetaSigmaDerivatives(), 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}

References MathLib::v.

Referenced by ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianEvalConstitutiveSetting(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForPressureEquations(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector(), ProcessLib::TH2M::ConstitutiveRelations::PermeabilityModel< DisplacementDim >::eval(), ProcessLib::ThermoRichardsMechanics::PermeabilityModel< DisplacementDim >::eval(), fourthOrderRotationMatrix< 2 >(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getEpsilon(), ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::getEpsilon(), ProcessLib::getIntegrationPointKelvinVectorData(), MaterialLib::Solids::Ehlers::SolidEhlers< DisplacementDim >::getInternalVariables(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::Reflection::detail::GetFlattenedIPDataFromLocAsm< Dim, Accessor_IPDataVecInLocAsm, Accessor_CurrentLevelFromIPDataVecElement >::operator()(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::postTimestepConcreteWithBlockVectors(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::updateConstitutiveRelations(), and MaterialPropertyLib::updateVariableArrayValues().

◆ 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}

References MathLib::v.

◆ 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}

References MathLib::v.

◆ 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, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References kelvinVectorToTensor(), OGS_FATAL, and MathLib::v.

◆ 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.

◆ KMnan()

template<int DisplacementDim>
auto MathLib::KelvinVector::KMnan ( )
constexpr

Returns an expressions for a Kelvin matrix filled with NaN.

Definition at line 83 of file KelvinVector.h.

84{
86 std::numeric_limits<double>::quiet_NaN());
87}

◆ KMzero()

template<int DisplacementDim>
auto MathLib::KelvinVector::KMzero ( )
constexpr

Returns an expressions for a Kelvin matrix filled with zero.

Definition at line 68 of file KelvinVector.h.

◆ KVnan()

template<int DisplacementDim>
auto MathLib::KelvinVector::KVnan ( )
constexpr

Returns an expressions for a Kelvin vector filled with NaN.

Definition at line 75 of file KelvinVector.h.

76{
78 std::numeric_limits<double>::quiet_NaN());
79}
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType

◆ KVzero()

template<int DisplacementDim>
auto MathLib::KelvinVector::KVzero ( )
constexpr

Returns an expressions for a Kelvin vector filled with zero.

Definition at line 61 of file KelvinVector.h.

◆ liftVectorToKelvin()

template<int DisplacementDim>
Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> MathLib::KelvinVector::liftVectorToKelvin ( Eigen::Matrix< double, DisplacementDim, 1 > const & v)

Lifting of a vector to a Kelvin matrix. \( a -> A\) s.t. \(k_{ij} a_{j} = A_[\alpha\beta} k_{\beta} \). Conversion for 2D -> 4D and 3D -> 6D are implemented.

Referenced by ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian().

◆ liftVectorToKelvin< 2 >()

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

Definition at line 162 of file KelvinVector.cpp.

186{
187 Eigen::Matrix<double, 2, 4> m;
188 m << v[0], 0, 0, v[1] / std::sqrt(2.), 0, v[1], 0, v[0] / std::sqrt(2.);
189 return m;
190}

◆ liftVectorToKelvin< 3 >()

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

Definition at line 162 of file KelvinVector.cpp.

195{
196 Eigen::Matrix<double, 3, 6> m;
197 m << v[0], 0, 0, v[1] / std::sqrt(2.), 0, v[2] / std::sqrt(2.), 0, v[1], 0,
198 v[0] / std::sqrt(2.), v[2] / std::sqrt(2.), 0, 0, 0, v[2], 0,
199 v[1] / std::sqrt(2.), v[0] / std::sqrt(2.);
200 return m;
201}

◆ reduceKelvinToVector()

template<int DisplacementDim>
Eigen::Matrix< double, DisplacementDim, 1 > MathLib::KelvinVector::reduceKelvinToVector ( Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> const & m)

Reducing a Kelvin matrix to a vector. Conversion for 4D -> 2D and 6D -> 3D are implemented.

◆ reduceKelvinToVector< 2 >()

template<>
Eigen::Matrix< double, 2, 1 > MathLib::KelvinVector::reduceKelvinToVector< 2 > ( Eigen::Matrix< double, 2, 4 > const & m)

Definition at line 162 of file KelvinVector.cpp.

206{
207 assert(m(0, 1) == 0);
208 assert(m(0, 2) == 0);
209 assert(m(1, 0) == 0);
210 assert(m(1, 2) == 0);
211 Eigen::Matrix<double, 2, 1> v;
212 v << m(0, 0), m(1, 1);
213 return v;
214}

◆ reduceKelvinToVector< 3 >()

template<>
Eigen::Matrix< double, 3, 1 > MathLib::KelvinVector::reduceKelvinToVector< 3 > ( Eigen::Matrix< double, 3, 6 > const & m)

Definition at line 162 of file KelvinVector.cpp.

219{
220 assert(m(0, 1) == 0);
221 assert(m(0, 2) == 0);
222 assert(m(0, 4) == 0);
223 assert(m(1, 0) == 0);
224 assert(m(1, 2) == 0);
225 assert(m(1, 5) == 0);
226 assert(m(2, 0) == 0);
227 assert(m(2, 1) == 0);
228 assert(m(2, 3) == 0);
229 assert(std::abs(m(0, 3) - m(2, 4)) <
230 std::numeric_limits<double>::epsilon());
231 assert(std::abs(m(0, 5) - m(1, 4)) <
232 std::numeric_limits<double>::epsilon());
233 assert(std::abs(m(1, 3) - m(2, 5)) <
234 std::numeric_limits<double>::epsilon());
235 Eigen::Matrix<double, 3, 1> v;
236 v << m(0, 0), m(1, 1), m(2, 2);
237 return v;
238}

◆ 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 201 of file KelvinVector.h.

202{
203 static_assert(
204 (Eigen::MatrixBase<Derived>::ColsAtCompileTime == 1) ||
205 (Eigen::MatrixBase<Derived>::ColsAtCompileTime == Eigen::Dynamic),
206 "KelvinVector must be a column vector");
207 if (v.cols() != 1)
208 {
209 OGS_FATAL(
210 "KelvinVector must be a column vector, but input has {:d} columns.",
211 v.cols());
212 }
213
214 Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1>
215 result;
216 if (v.rows() == 4)
217 {
218 result.resize(4, 1);
219 result << v[0], v[1], v[2], v[3] * std::sqrt(2.);
220 }
221 else if (v.rows() == 6)
222 {
223 result.resize(6, 1);
224 result << v[0], v[1], v[2], v[3] * std::sqrt(2.), v[4] * std::sqrt(2.),
225 v[5] * std::sqrt(2.);
226 }
227 else
228 {
229 OGS_FATAL(
230 "Symmetric tensor to Kelvin vector conversion expected an input "
231 "vector of size 4 or 6, but a vector of size {:d} was given.",
232 v.size());
233 }
234 return result;
235}

References OGS_FATAL, and MathLib::v.

Referenced by ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::PhaseField::PhaseFieldLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::TH2M::TH2MLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::initializeConcrete(), ProcessLib::setIntegrationPointKelvinVectorData(), ProcessLib::Reflection::detail::setIPData(), 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 242 of file KelvinVector.h.

244{
245 constexpr int kelvin_vector_size =
246 kelvin_vector_dimensions(DisplacementDim);
247
248 if (values.size() != kelvin_vector_size)
249 {
250 OGS_FATAL(
251 "Symmetric tensor to Kelvin vector conversion expected an input "
252 "vector of size {:d}, but a vector of size {:d} was given.",
253 kelvin_vector_size, values.size());
254 }
255
256 return symmetricTensorToKelvinVector(
257 Eigen::Map<typename MathLib::KelvinVector::KelvinVectorType<
258 DisplacementDim> const>(
259 values.data(), kelvin_vector_dimensions(DisplacementDim), 1));
260}

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

◆ tensorToKelvin()

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

◆ 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
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
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}