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
111 static Eigen::Matrix<double, KelvinVectorSize, 1> const ones2;
112
114 static double determinant(
115 Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
116
119 static double equivalentStress(
120 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
121
123 static double FrobeniusNorm(
124 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
125
128 static double J2(
129 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
130
133 static double J3(
134 Eigen::Matrix<double, KelvinVectorSize, 1> const& deviatoric_v);
135
137 static double trace(Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
138
141 static Eigen::Vector3d diagonal(
142 Eigen::Matrix<double, KelvinVectorSize, 1> const& v);
143};
144
145//
146// Inverses of a Kelvin vector.
147//
148
151template <int KelvinVectorSize>
152Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1>
153inverse(Eigen::Matrix<double,
154 KelvinVectorSize,
155 1,
156 Eigen::ColMajor,
157 KelvinVectorSize,
158 1> const& v);
159
162template <int KelvinVectorSize>
163Eigen::Matrix<double, 3, 3> kelvinVectorToTensor(Eigen::Matrix<double,
164 KelvinVectorSize,
165 1,
166 Eigen::ColMajor,
167 KelvinVectorSize,
168 1> const& v);
169
172template <int DisplacementDim>
174 Eigen::Matrix<double, 3, 3> const& m);
175
187template <int KelvinVectorSize>
188Eigen::Matrix<double, KelvinVectorSize, 1, Eigen::ColMajor, KelvinVectorSize, 1>
189kelvinVectorToSymmetricTensor(Eigen::Matrix<double,
190 KelvinVectorSize,
191 1,
192 Eigen::ColMajor,
193 KelvinVectorSize,
194 1> const& v);
195
203template <typename Derived>
204Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1>
205symmetricTensorToKelvinVector(Eigen::MatrixBase<Derived> const& v)
206{
207 static_assert(
208 (Eigen::MatrixBase<Derived>::ColsAtCompileTime == 1) ||
209 (Eigen::MatrixBase<Derived>::ColsAtCompileTime == Eigen::Dynamic),
210 "KelvinVector must be a column vector");
211 if (v.cols() != 1)
212 {
213 OGS_FATAL(
214 "KelvinVector must be a column vector, but input has {:d} columns.",
215 v.cols());
216 }
217
218 Eigen::Matrix<double, Eigen::MatrixBase<Derived>::RowsAtCompileTime, 1>
219 result;
220 if (v.rows() == 4)
221 {
222 result.resize(4, 1);
223 result << v[0], v[1], v[2], v[3] * std::sqrt(2.);
224 }
225 else if (v.rows() == 6)
226 {
227 result.resize(6, 1);
228 result << v[0], v[1], v[2], v[3] * std::sqrt(2.), v[4] * std::sqrt(2.),
229 v[5] * std::sqrt(2.);
230 }
231 else
232 {
233 OGS_FATAL(
234 "Symmetric tensor to Kelvin vector conversion expected an input "
235 "vector of size 4 or 6, but a vector of size {:d} was given.",
236 v.size());
237 }
238 return result;
239}
240
245template <int DisplacementDim>
247 std::vector<double> const& values)
248{
249 constexpr int kelvin_vector_size =
250 kelvin_vector_dimensions(DisplacementDim);
251
252 if (values.size() != kelvin_vector_size)
253 {
254 OGS_FATAL(
255 "Symmetric tensor to Kelvin vector conversion expected an input "
256 "vector of size {:d}, but a vector of size {:d} was given.",
257 kelvin_vector_size, values.size());
258 }
259
261 Eigen::Map<typename MathLib::KelvinVector::KelvinVectorType<
262 DisplacementDim> const>(
263 values.data(), kelvin_vector_dimensions(DisplacementDim), 1));
264}
265
269template <int DisplacementDim>
270Eigen::Matrix<double, DisplacementDim,
271 kelvin_vector_dimensions(DisplacementDim)>
272liftVectorToKelvin(Eigen::Matrix<double, DisplacementDim, 1> const& v);
273
276template <int DisplacementDim>
277Eigen::Matrix<double, DisplacementDim, 1> reduceKelvinToVector(
278 Eigen::Matrix<double, DisplacementDim,
279 kelvin_vector_dimensions(DisplacementDim)> const& m);
280
285template <int DisplacementDim>
287 Eigen::Matrix<double, DisplacementDim, DisplacementDim, Eigen::ColMajor,
288 DisplacementDim, DisplacementDim> const& transformation);
289
290} // namespace KelvinVector
291} // namespace MathLib
292
293#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 ones2
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)