Newer
Older
* Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#include "KelvinVector.h"
#include "BaseLib/Error.h"
namespace KelvinVector
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
{
template <>
double Invariants<6>::determinant(Eigen::Matrix<double, 6, 1> const& v)
{
return v(0) * v(1) * v(2) + v(3) * v(4) * v(5) / std::sqrt(2.) -
v(3) * v(3) * v(2) / 2. - v(4) * v(4) * v(0) / 2. -
v(5) * v(5) * v(1) / 2.;
}
template <>
double Invariants<4>::determinant(Eigen::Matrix<double, 4, 1> const& v)
{
return v(2) * (v(0) * v(1) - v(3) * v(3) / 2.);
}
template <>
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> inverse(
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> const& v)
{
assert(Invariants<4>::determinant(v) != 0);
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> inv;
inv(0) = v(1) * v(2);
inv(1) = v(0) * v(2);
inv(2) = v(0) * v(1) - v(3) * v(3) / 2.;
inv(3) = -v(3) * v(2);
return inv / Invariants<4>::determinant(v);
}
template <>
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> inverse(
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> const& v)
{
assert(Invariants<6>::determinant(v) != 0);
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> inv;
inv(0) = v(1) * v(2) - v(4) * v(4) / 2.;
inv(1) = v(0) * v(2) - v(5) * v(5) / 2.;
inv(2) = v(0) * v(1) - v(3) * v(3) / 2.;
inv(3) = v(4) * v(5) / std::sqrt(2.) - v(3) * v(2);
inv(4) = v(3) * v(5) / std::sqrt(2.) - v(4) * v(0);
inv(5) = v(4) * v(3) / std::sqrt(2.) - v(1) * v(5);
return inv / Invariants<6>::determinant(v);
}
Eigen::Matrix<double, 3, 3> kelvinVectorToTensor(
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> const& v)
{
Eigen::Matrix<double, 3, 3> m;
m << v[0], v[3] / std::sqrt(2.), 0, v[3] / std::sqrt(2.), v[1], 0, 0, 0,
v[2];
return m;
}
template <>
Eigen::Matrix<double, 3, 3> kelvinVectorToTensor(
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> const& v)
{
Eigen::Matrix<double, 3, 3> m;
m << v[0], v[3] / std::sqrt(2.), v[5] / std::sqrt(2.), v[3] / std::sqrt(2.),
v[1], v[4] / std::sqrt(2.), v[5] / std::sqrt(2.), v[4] / std::sqrt(2.),
v[2];
return m;
}
template <>
Eigen::Matrix<double, 3, 3> kelvinVectorToTensor(Eigen::Matrix<double,
Eigen::Dynamic,
1,
Eigen::ColMajor,
Eigen::Dynamic,
1> const& v)
{
if (v.size() == 4)
{
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> v4;
v4 << v[0], v[1], v[2], v[3];
return kelvinVectorToTensor(v4);
}
if (v.size() == 6)
{
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> v6;
v6 << v[0], v[1], v[2], v[3], v[4], v[5];
return kelvinVectorToTensor(v6);
}
OGS_FATAL(
"Conversion of dynamic Kelvin vector of size {:d} to a tensor is not "
"possible. Kelvin vector must be of size 4 or 6.",
v.size());
}
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
template <>
KelvinVectorType<2> tensorToKelvin<2>(Eigen::Matrix<double, 3, 3> const& m)
{
assert(std::abs(m(0, 1) - m(1, 0)) <
std::numeric_limits<double>::epsilon());
assert(m(0, 2) == 0);
assert(m(1, 2) == 0);
assert(m(2, 0) == 0);
assert(m(2, 1) == 0);
KelvinVectorType<2> v;
v << m(0, 0), m(1, 1), m(2, 2), m(0, 1) * std::sqrt(2.);
return v;
}
template <>
KelvinVectorType<3> tensorToKelvin<3>(Eigen::Matrix<double, 3, 3> const& m)
{
assert(std::abs(m(0, 1) - m(1, 0)) <
std::numeric_limits<double>::epsilon());
assert(std::abs(m(1, 2) - m(2, 1)) <
std::numeric_limits<double>::epsilon());
assert(std::abs(m(0, 2) - m(2, 0)) <
std::numeric_limits<double>::epsilon());
KelvinVectorType<3> v;
v << m(0, 0), m(1, 1), m(2, 2), m(0, 1) * std::sqrt(2.),
m(1, 2) * std::sqrt(2.), m(0, 2) * std::sqrt(2.);
return v;
}
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
template <>
Eigen::Matrix<double, 4, 1> kelvinVectorToSymmetricTensor(
Eigen::Matrix<double, 4, 1, Eigen::ColMajor, 4, 1> const& v)
{
Eigen::Matrix<double, 4, 1> m;
m << v[0], v[1], v[2], v[3] / std::sqrt(2.);
return m;
}
template <>
Eigen::Matrix<double, 6, 1> kelvinVectorToSymmetricTensor(
Eigen::Matrix<double, 6, 1, Eigen::ColMajor, 6, 1> const& v)
{
Eigen::Matrix<double, 6, 1> m;
m << v[0], v[1], v[2], v[3] / std::sqrt(2.), v[4] / std::sqrt(2.),
v[5] / std::sqrt(2.);
return m;
}
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)
{
if (v.size() == 4)
{
return kelvinVectorToSymmetricTensor<4>(v);
}
{
return kelvinVectorToSymmetricTensor<6>(v);
}
OGS_FATAL(
"Kelvin vector to tensor conversion expected an input vector of size 4 "
template <>
KelvinMatrixType<2> fourthOrderRotationMatrix<2>(
Eigen::Matrix<double, 2, 2, Eigen::ColMajor, 2, 2> const& transformation)
{
// 1-based index access for convenience.
auto Q = [&](int const i, int const j) {
return transformation(i - 1, j - 1);
};
MathLib::KelvinVector::KelvinMatrixType<2> R;
R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), 0,
std::sqrt(2) * Q(1, 1) * Q(1, 2), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
0, std::sqrt(2) * Q(2, 1) * Q(2, 2), 0, 0, 1, 0,
std::sqrt(2) * Q(1, 1) * Q(2, 1), std::sqrt(2) * Q(1, 2) * Q(2, 2), 0,
Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1);
return R;
}
template <>
KelvinMatrixType<3> fourthOrderRotationMatrix<3>(
Eigen::Matrix<double, 3, 3, Eigen::ColMajor, 3, 3> const& transformation)
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
{
// 1-based index access for convenience.
auto Q = [&](int const i, int const j) {
return transformation(i - 1, j - 1);
};
MathLib::KelvinVector::KelvinMatrixType<3> R;
R << Q(1, 1) * Q(1, 1), Q(1, 2) * Q(1, 2), Q(1, 3) * Q(1, 3),
std::sqrt(2) * Q(1, 1) * Q(1, 2), std::sqrt(2) * Q(1, 2) * Q(1, 3),
std::sqrt(2) * Q(1, 1) * Q(1, 3), Q(2, 1) * Q(2, 1), Q(2, 2) * Q(2, 2),
Q(2, 3) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 2),
std::sqrt(2) * Q(2, 2) * Q(2, 3), std::sqrt(2) * Q(2, 1) * Q(2, 3),
Q(3, 1) * Q(3, 1), Q(3, 2) * Q(3, 2), Q(3, 3) * Q(3, 3),
std::sqrt(2) * Q(3, 1) * Q(3, 2), std::sqrt(2) * Q(3, 2) * Q(3, 3),
std::sqrt(2) * Q(3, 1) * Q(3, 3), std::sqrt(2) * Q(1, 1) * Q(2, 1),
std::sqrt(2) * Q(1, 2) * Q(2, 2), std::sqrt(2) * Q(1, 3) * Q(2, 3),
Q(1, 1) * Q(2, 2) + Q(1, 2) * Q(2, 1),
Q(1, 2) * Q(2, 3) + Q(1, 3) * Q(2, 2),
Q(1, 1) * Q(2, 3) + Q(1, 3) * Q(2, 1), std::sqrt(2) * Q(2, 1) * Q(3, 1),
std::sqrt(2) * Q(2, 2) * Q(3, 2), std::sqrt(2) * Q(2, 3) * Q(3, 3),
Q(2, 1) * Q(3, 2) + Q(2, 2) * Q(3, 1),
Q(2, 2) * Q(3, 3) + Q(2, 3) * Q(3, 2),
Q(2, 1) * Q(3, 3) + Q(2, 3) * Q(3, 1), std::sqrt(2) * Q(1, 1) * Q(3, 1),
std::sqrt(2) * Q(1, 2) * Q(3, 2), std::sqrt(2) * Q(1, 3) * Q(3, 3),
Q(1, 1) * Q(3, 2) + Q(1, 2) * Q(3, 1),
Q(1, 2) * Q(3, 3) + Q(1, 3) * Q(3, 2),
Q(1, 1) * Q(3, 3) + Q(1, 3) * Q(3, 1);
return R;
}
} // namespace KelvinVector
} // namespace MathLib