OGS
KelvinVector.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <Eigen/Core>
13
14#include "BaseLib/Error.h"
15
16namespace MathLib
17{
21namespace KelvinVector
22{
24constexpr int kelvin_vector_dimensions(int const displacement_dim)
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}
38
39//
40// Kelvin vector and matrix templates for given displacement dimension.
41//
42
46template <int DisplacementDim>
48 Eigen::Matrix<double, kelvin_vector_dimensions(DisplacementDim), 1,
49 Eigen::ColMajor>;
50
54template <int DisplacementDim>
56 Eigen::Matrix<double, kelvin_vector_dimensions(DisplacementDim),
57 kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor>;
58
60template <int DisplacementDim>
61constexpr auto KVzero()
62{
64}
65
67template <int DisplacementDim>
68constexpr auto KMzero()
69{
71}
72
74template <int DisplacementDim>
75constexpr auto KVnan()
76{
78 std::numeric_limits<double>::quiet_NaN());
79}
80
82template <int DisplacementDim>
83constexpr auto KMnan()
84{
86 std::numeric_limits<double>::quiet_NaN());
87}
88
92template <int KelvinVectorSize>
93struct Invariants final
94{
95 static_assert(KelvinVectorSize == 4 || KelvinVectorSize == 6,
96 "KelvinVector invariants for vectors of size different than "
97 "4 or 6 is not allowed.");
100 static Eigen::Matrix<double, KelvinVectorSize, KelvinVectorSize> const
104 static Eigen::Matrix<double, KelvinVectorSize, KelvinVectorSize> const
107 static Eigen::Matrix<double, KelvinVectorSize, 1> const identity2;
108
110 static double determinant(
111 Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
112
115 static double equivalentStress(
116 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
117
119 static double FrobeniusNorm(
120 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
121
124 static double J2(
125 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
126
129 static double J3(
130 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
131
133 static double trace(Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
134
137 static Eigen::Vector3d diagonal(
138 Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
139};
140
141//
142// Inverses of a Kelvin vector.
143//
144
147template <int KelvinVectorSize>
148Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1>
149inverse(Eigen::Matrix<double,
150 KelvinVectorSize,
151 1,
152 Eigen::ColMajor,
153 KelvinVectorSize,
154 1> const& v);
155
158template <int KelvinVectorSize>
159Eigen::Matrix<double, 3, 3> kelvinVectorToTensor(Eigen::Matrix<double,
160 KelvinVectorSize,
161 1,
162 Eigen::ColMajor,
163 KelvinVectorSize,
164 1> const& v);
165
168template <int DisplacementDim>
170 Eigen::Matrix<double, 3, 3> const& m);
171
183template <int KelvinVectorSize>
184Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1>
185kelvinVectorToSymmetricTensor(Eigen::Matrix<double,
186 KelvinVectorSize,
187 1,
188 Eigen::ColMajor,
189 KelvinVectorSize,
190 1> const& v);
191
199template <typename Derived>
200Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1>
201symmetricTensorToKelvinVector(Eigen::MatrixBase<Derived> const& v)
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}
236
241template <int DisplacementDim>
243 std::vector<double> const& values)
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
257 Eigen::Map<typename MathLib::KelvinVector::KelvinVectorType<
258 DisplacementDim> const>(
259 values.data(), kelvin_vector_dimensions(DisplacementDim), 1));
260}
261
265template <int DisplacementDim>
266Eigen::Matrix<double, DisplacementDim,
267 kelvin_vector_dimensions(DisplacementDim)>
268liftVectorToKelvin(Eigen::Matrix<double, DisplacementDim, 1> const& v);
269
272template <int DisplacementDim>
273Eigen::Matrix<double, DisplacementDim, 1> reduceKelvinToVector(
274 Eigen::Matrix<double, DisplacementDim,
275 kelvin_vector_dimensions(DisplacementDim)> const& m);
276
281template <int DisplacementDim>
283 Eigen::Matrix<double, DisplacementDim, DisplacementDim, Eigen::ColMajor,
284 DisplacementDim, DisplacementDim> const& transformation);
285
286} // namespace KelvinVector
287} // namespace MathLib
288
289#include "KelvinVector-impl.h"
#define OGS_FATAL(...)
Definition Error.h:26
Eigen::Matrix< double, DisplacementDim, 1 > reduceKelvinToVector(Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> const &m)
KelvinMatrixType< DisplacementDim > fourthOrderRotationMatrix(Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::ColMajor, DisplacementDim, DisplacementDim > const &transformation)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
constexpr auto KVnan()
Returns an expressions for a Kelvin vector filled with NaN.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr auto KMzero()
Returns an expressions for a Kelvin matrix filled with zero.
constexpr auto KVzero()
Returns an expressions for a Kelvin vector filled with zero.
Eigen::Matrix< double, DisplacementDim, kelvin_vector_dimensions(DisplacementDim)> liftVectorToKelvin(Eigen::Matrix< double, DisplacementDim, 1 > const &v)
constexpr auto KMnan()
Returns an expressions for a Kelvin matrix filled with NaN.
KelvinVectorType< DisplacementDim > tensorToKelvin(Eigen::Matrix< double, 3, 3 > const &m)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > inverse(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
static const double v
KV::KelvinVectorType< DisplacementDim > KelvinVector
Definition Base.h:27
static double FrobeniusNorm(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
Get the norm of the deviatoric stress.
static Eigen::Matrix< double, KelvinVectorSize, 1 > const identity2
Kelvin mapping of 2nd order identity tensor.
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 Eigen::Matrix< double, KelvinVectorSize, KelvinVectorSize > const spherical_projection
static double J2(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)
static Eigen::Matrix< double, KelvinVectorSize, KelvinVectorSize > const deviatoric_projection
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 J3(Eigen::Matrix< double, KelvinVectorSize, 1 > const &deviatoric_v)