version 0.4.1
tensorutils.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10
11#include <numeric>
12#include <ranges>
13#include <unsupported/Eigen/CXX11/Tensor>
14
15#include <dune/common/promotiontraits.hh>
16
18#include <ikarus/utils/math.hh>
19namespace Ikarus {
20
31template <typename Derived, typename T, auto rank>
32Eigen::Tensor<typename Derived::Scalar, rank> tensorView(const Eigen::EigenBase<Derived>& matrix,
33 const std::array<T, rank>& dims) {
34 return Eigen::TensorMap<const Eigen::TensorFixedSize<
35 const typename Derived::Scalar, Eigen::Sizes<Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>>>(
36 matrix.derived().eval().data(), dims);
37}
38
47auto dyadic(const auto& A_ij, const auto& B_kl) {
48 Eigen::array<Eigen::IndexPair<long>, 0> empty_index_list = {};
49 return A_ij.contract(B_kl, empty_index_list).eval();
50}
51
60template <typename ST, int size>
61auto dyadic(const Eigen::Vector<ST, size>& a, const Eigen::Vector<ST, size>& b) {
62 return (a * b.transpose()).eval();
63}
64
72template <typename ScalarType = double, int dim = 3>
74 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> idTensor;
75 for (int i = 0; i < dim; ++i)
76 for (int j = 0; j < dim; ++j)
77 for (int k = 0; k < dim; ++k)
78 for (int l = 0; l < dim; ++l)
79 idTensor(i, j, k, l) = 0.5 * ((i == k) * (j == l) + (i == l) * (j == k));
80 return idTensor;
81}
82
94template <typename ScalarType = double, int dim = 3>
95auto symmetricFourthOrder(const auto& A, const auto& B) {
96 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> res;
97 for (int i = 0; i < dim; ++i)
98 for (int j = 0; j < dim; ++j)
99 for (int k = 0; k < dim; ++k)
100 for (int l = 0; l < dim; ++l)
101 res(i, j, k, l) = 0.5 * (A(i, k) * B(j, l) + A(i, l) * B(j, k));
102 return res;
103}
104
113template <typename ScalarType = double, int dim = 3>
115 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> idTensor;
116 idTensor.setZero();
117 for (int i = 0; i < dim; ++i)
118 for (int k = 0; k < dim; ++k)
119 idTensor(i, i, k, k) = 1.0;
120 return idTensor;
121}
122
134template <typename AType, typename BType>
135auto fourthOrderIKJL(const Eigen::MatrixBase<AType>& A, const Eigen::MatrixBase<BType>& B) {
136 static_assert(AType::RowsAtCompileTime == BType::RowsAtCompileTime);
137 static_assert(AType::ColsAtCompileTime == BType::ColsAtCompileTime);
138 using ScalarResultType = typename Dune::PromotionTraits<typename AType::Scalar, typename BType::Scalar>::PromotedType;
139 constexpr int dim = AType::RowsAtCompileTime;
140 Eigen::TensorFixedSize<ScalarResultType, Eigen::Sizes<dim, dim, dim, dim>> res;
141 for (int i = 0; i < dim; ++i)
142 for (int j = 0; j < dim; ++j)
143 for (int k = 0; k < dim; ++k)
144 for (int l = 0; l < dim; ++l)
145 res(i, j, k, l) = A(i, k) * B(j, l);
146 return res;
147}
148
157template <typename ScalarType, long int dim>
158auto symTwoSlots(const Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>>& t,
159 const std::array<size_t, 2>& slots) {
160 std::array<size_t, 4> shuffleSlot;
161 std::iota(shuffleSlot.begin(), shuffleSlot.end(), 0); // create 0,1,2,3 array
162 std::swap(shuffleSlot[slots[0]], shuffleSlot[slots[1]]); // swap the given slots
163 return (0.5 * (t + t.shuffle(shuffleSlot))).eval();
164}
165
180template <int dim = 3>
181requires(dim == 2 or dim == 3)
182constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept {
183 if constexpr (dim == 2) {
184 if (i == j) // _00 -> 0, _11 -> 1
185 return i;
186 if ((i == 0 and j == 1) or (i == 1 and j == 0)) // _01 and _10 --> 2
187 return 2;
188 assert(i < 2 and j < 2 && "For Voigt notation the indices need to be 0 or 1.");
189 __builtin_unreachable();
190 } else {
191 if (i == j) // _00 -> 0, _11 -> 1, _22 -> 2
192 return i;
193 if ((i == 1 and j == 2) or (i == 2 and j == 1)) // _12 and _21 --> 3
194 return 3;
195 if ((i == 0 and j == 2) or (i == 2 and j == 0)) // _02 and _20 --> 4
196 return 4;
197 if ((i == 0 and j == 1) or (i == 1 and j == 0)) // _01 and _10 --> 5
198 return 5;
199 assert(i < 3 and j < 3 && "For Voigt notation the indices need to be 0,1 or 2.");
200 __builtin_unreachable();
201 }
202}
203
218template <typename ScalarType = double>
219Eigen::Matrix<ScalarType, 6, 6> toVoigt(const Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3, 3, 3>>& ft) {
220 Eigen::Matrix<ScalarType, 6, 6> mat;
221 for (Eigen::Index i = 0; i < 3; ++i)
222 for (Eigen::Index j = 0; j < 3; ++j)
223 for (Eigen::Index k = 0; k < 3; ++k)
224 for (Eigen::Index l = 0; l < 3; ++l)
225 mat(toVoigt(i, j), toVoigt(k, l)) = ft(i, j, k, l);
226 return mat;
227}
228
248template <typename ST, int size, int Options, int maxSize>
249requires((size > 0 and size <= 3) or (maxSize > 0 and maxSize <= 3 and size == Eigen::Dynamic))
250auto toVoigt(const Eigen::Matrix<ST, size, size, Options, maxSize, maxSize>& E, bool isStrain = true) {
251 constexpr bool isFixedSized = (size != Eigen::Dynamic);
252 const ST possibleStrainFactor = isStrain ? 2.0 : 1.0;
253
254 const size_t inputSize = isFixedSized ? size : E.rows();
255 auto EVoigt = [&]() {
256 if constexpr (isFixedSized) {
257 Eigen::Vector<ST, (size * (size + 1)) / 2> EVoigt;
258 EVoigt.template head<size>() = E.diagonal();
259 return EVoigt;
260 } else {
261 Eigen::Matrix<ST, Eigen::Dynamic, 1, Options, 6, 1> EVoigt;
262 EVoigt.resize((inputSize * (inputSize + 1)) / 2);
263 EVoigt.template head(inputSize) = E.diagonal();
264 return EVoigt;
265 }
266 }();
267
268 if (inputSize == 2)
269 EVoigt(2) = E(0, 1) * possibleStrainFactor;
270 else if (inputSize == 3) {
271 EVoigt(inputSize) = E(1, 2) * possibleStrainFactor;
272 EVoigt(inputSize + 1) = E(0, 2) * possibleStrainFactor;
273 EVoigt(inputSize + 2) = E(0, 1) * possibleStrainFactor;
274 }
275 return EVoigt;
276}
277
293template <typename ST, int size, int Options, int maxSize>
294requires((size == 1 or size == 3 or size == 6) or
295 ((maxSize == 1 or maxSize == 3 or maxSize == 6) and size == Eigen::Dynamic))
296auto fromVoigt(const Eigen::Matrix<ST, size, 1, Options, maxSize, 1>& EVoigt, bool isStrain = true) {
297 constexpr bool isFixedSized = (size != Eigen::Dynamic);
298 const ST possibleStrainFactor = isStrain ? 0.5 : 1.0;
299
300 const size_t inputSize = isFixedSized ? size : EVoigt.size();
301 const size_t matrixSize =
302 isFixedSized ? (-1 + ct_sqrt(1 + 8 * size)) / 2 : (-1 + static_cast<int>(std::sqrt(1 + 8 * inputSize))) / 2;
303
304 auto E = [&]() {
305 if constexpr (isFixedSized) {
306 Eigen::Matrix<ST, matrixSize, matrixSize> E;
307 E.diagonal() = EVoigt.template head<matrixSize>();
308 return E;
309 } else {
310 Eigen::Matrix<ST, Eigen::Dynamic, Eigen::Dynamic, Options, 3, 3> E;
311 E.resize(matrixSize, matrixSize);
312 E.diagonal() = EVoigt.template head(matrixSize);
313 return E;
314 }
315 }();
316
317 if (matrixSize == 2) {
318 E(0, 1) = E(1, 0) = EVoigt(2) * possibleStrainFactor;
319 } else if (matrixSize == 3) {
320 E(2, 1) = E(1, 2) = EVoigt(matrixSize) * possibleStrainFactor;
321 E(2, 0) = E(0, 2) = EVoigt(matrixSize + 1) * possibleStrainFactor;
322 E(1, 0) = E(0, 1) = EVoigt(matrixSize + 2) * possibleStrainFactor;
323 }
324
325 return E;
326}
327
341constexpr std::array<size_t, 2> fromVoigt(size_t i) {
342 if (i < 3) // _00 -> 0, _11 -> 1, _22 -> 2
343 return {i, i};
344 else if (i == 3)
345 return {1, 2};
346 else if (i == 4)
347 return {0, 2};
348 else if (i == 5)
349 return {0, 1};
350 else {
351 assert(i < 6 && "For Voigt notation the indices need to be 0 and 5.");
352 __builtin_unreachable();
353 }
354}
355
367template <typename ScalarType>
368auto fromVoigt(const Eigen::Matrix<ScalarType, 6, 6>& CVoigt) {
369 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3, 3, 3>> C;
370 // size_t iR=0,jR=0;
371 for (size_t i = 0; i < 6; ++i) {
372 for (size_t j = 0; j < 6; ++j) {
373 auto firstIndices = fromVoigt(i);
374 auto secondIndices = fromVoigt(j);
375 C(firstIndices[0], firstIndices[1], secondIndices[0], secondIndices[1]) = CVoigt(i, j);
376 }
377 }
378 return C;
379}
380
392template <typename GEO>
393requires(GEO::mydimension == 2)
394Eigen::Matrix3d transformationMatrix(const GEO& geometry, const Dune::FieldVector<double, 2>& pos) {
395 const auto jacobian = toEigen(geometry.jacobianTransposed(pos)).eval();
396
397 const auto J11 = jacobian(0, 0);
398 const auto J12 = jacobian(0, 1);
399 const auto J21 = jacobian(1, 0);
400 const auto J22 = jacobian(1, 1);
401
402 Eigen::Matrix3d T{
403 { J11 * J11, J12 * J12, J11 * J12},
404 { J21 * J21, J22 * J22, J21 * J22},
405 {2.0 * J11 * J21, 2.0 * J12 * J22, J21 * J12 + J11 * J22}
406 };
407
408 return T;
409}
410
422template <typename GEO>
423requires(GEO::mydimension == 3)
424Eigen::Matrix<double, 6, 6> transformationMatrix(const GEO& geometry, const Dune::FieldVector<double, 3>& pos) {
425 const auto jacobian = toEigen(geometry.jacobianTransposed(pos)).eval();
426
427 const auto J11 = jacobian(0, 0);
428 const auto J12 = jacobian(0, 1);
429 const auto J13 = jacobian(0, 2);
430 const auto J21 = jacobian(1, 0);
431 const auto J22 = jacobian(1, 1);
432 const auto J23 = jacobian(1, 2);
433 const auto J31 = jacobian(2, 0);
434 const auto J32 = jacobian(2, 1);
435 const auto J33 = jacobian(2, 2);
436
437 // clang-format off
438 Eigen::Matrix<double, 6, 6> T {
439 {J11 * J11, J12 * J12, J13 * J13, J11 * J12, J11 * J13, J12 * J13},
440 {J21 * J21, J22 * J22, J23 * J23, J21 * J22, J21 * J23, J22 * J23},
441 {J31 * J31, J32 * J32, J33 * J33, J31 * J32, J31 * J33, J32 * J33},
442 {2.0 * J11 * J21, 2.0 * J12 * J22, 2.0 * J13 * J23, J11 * J22 + J21 * J12, J11 * J23 + J21 * J13, J12 * J23 + J22 * J13},
443 {2.0 * J11 * J31, 2.0 * J12 * J32, 2.0 * J13 * J33, J11 * J32 + J31 * J12, J11 * J33 + J31 * J13, J12 * J33 + J32 * J13},
444 {2.0 * J31 * J21, 2.0 * J32 * J22, 2.0 * J33 * J23, J31 * J22 + J21 * J32, J31 * J23 + J21 * J33, J32 * J23 + J22 * J33}
445 };
446 // clang-format on
447
448 return T;
449}
450
451} // namespace Ikarus
Implementation of math related algorithms.
auto symmetricIdentityFourthOrder()
Generates a symmetric identity fourth-order tensor.
Definition: tensorutils.hh:73
auto symmetricFourthOrder(const auto &A, const auto &B)
Generates a symmetric fourth-order tensor based on two second-order input tensors.
Definition: tensorutils.hh:95
Eigen::Tensor< typename Derived::Scalar, rank > tensorView(const Eigen::EigenBase< Derived > &matrix, const std::array< T, rank > &dims)
View an Eigen matrix as an Eigen Tensor with specified dimensions.
Definition: tensorutils.hh:32
auto fourthOrderIKJL(const Eigen::MatrixBase< AType > &A, const Eigen::MatrixBase< BType > &B)
Computes the IKJL product of two matrices.
Definition: tensorutils.hh:135
auto dyadic(const auto &A_ij, const auto &B_kl)
Computes the dyadic product of two Eigen tensors.
Definition: tensorutils.hh:47
auto identityFourthOrder()
Generates an identity fourth-order tensor.
Definition: tensorutils.hh:114
auto fromVoigt(const Eigen::Matrix< ST, size, 1, Options, maxSize, 1 > &EVoigt, bool isStrain=true)
Converts a vector given in Voigt notation to a matrix.
Definition: tensorutils.hh:296
auto symTwoSlots(const Eigen::TensorFixedSize< ScalarType, Eigen::Sizes< dim, dim, dim, dim > > &t, const std::array< size_t, 2 > &slots)
Creates a symmetric fourth-order tensor in the two specified slots of the input tensor.
Definition: tensorutils.hh:158
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:182
Definition: assemblermanipulatorbuildingblocks.hh:22
Eigen::Matrix3d transformationMatrix(const GEO &geometry, const Dune::FieldVector< double, 2 > &pos)
Calculates the 2D transformation matrix.
Definition: tensorutils.hh:394
constexpr T ct_sqrt(T x)
Compile-time square root for integer types.
Definition: math.hh:47
Definition: truncatedconjugategradient.hh:24
Definition: utils/concepts.hh:30
Several concepts.