version 0.4
tensorutils.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 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
31 template <typename Derived, typename T, auto rank>
32 Eigen::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
47 auto 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
59 template <typename ScalarType = double, int dim = 3>
61 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> idTensor;
62 for (int i = 0; i < dim; ++i)
63 for (int j = 0; j < dim; ++j)
64 for (int k = 0; k < dim; ++k)
65 for (int l = 0; l < dim; ++l)
66 idTensor(i, j, k, l) = 0.5 * ((i == k) * (j == l) + (i == l) * (j == k));
67 return idTensor;
68 }
69
81 template <typename ScalarType = double, int dim = 3>
82 auto symmetricFourthOrder(const auto& A, const auto& B) {
83 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> res;
84 for (int i = 0; i < dim; ++i)
85 for (int j = 0; j < dim; ++j)
86 for (int k = 0; k < dim; ++k)
87 for (int l = 0; l < dim; ++l)
88 res(i, j, k, l) = 0.5 * (A(i, k) * B(j, l) + A(i, l) * B(j, k));
89 return res;
90 }
91
100 template <typename ScalarType = double, int dim = 3>
102 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>> idTensor;
103 idTensor.setZero();
104 for (int i = 0; i < dim; ++i)
105 for (int k = 0; k < dim; ++k)
106 idTensor(i, i, k, k) = 1.0;
107 return idTensor;
108 }
109
121 template <typename AType, typename BType>
122 auto fourthOrderIKJL(const Eigen::MatrixBase<AType>& A, const Eigen::MatrixBase<BType>& B) {
123 static_assert(AType::RowsAtCompileTime == BType::RowsAtCompileTime);
124 static_assert(AType::ColsAtCompileTime == BType::ColsAtCompileTime);
125 using ScalarResultType =
126 typename Dune::PromotionTraits<typename AType::Scalar, typename BType::Scalar>::PromotedType;
127 constexpr int dim = AType::RowsAtCompileTime;
128 Eigen::TensorFixedSize<ScalarResultType, Eigen::Sizes<dim, dim, dim, dim>> res;
129 for (int i = 0; i < dim; ++i)
130 for (int j = 0; j < dim; ++j)
131 for (int k = 0; k < dim; ++k)
132 for (int l = 0; l < dim; ++l)
133 res(i, j, k, l) = A(i, k) * B(j, l);
134 return res;
135 }
136
145 template <typename ScalarType, long int dim>
146 auto symTwoSlots(const Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<dim, dim, dim, dim>>& t,
147 const std::array<size_t, 2>& slots) {
148 std::array<size_t, 4> shuffleSlot;
149 std::iota(shuffleSlot.begin(), shuffleSlot.end(), 0); // create 0,1,2,3 array
150 std::swap(shuffleSlot[slots[0]], shuffleSlot[slots[1]]); // swap the given slots
151 return (0.5 * (t + t.shuffle(shuffleSlot))).eval();
152 }
153
167 constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept {
168 if (i == j) // _00 -> 0, _11 -> 1, _22 -> 2
169 return i;
170 if ((i == 1 and j == 2) or (i == 2 and j == 1)) // _12 and _21 --> 3
171 return 3;
172 if ((i == 0 and j == 2) or (i == 2 and j == 0)) // _02 and _20 --> 4
173 return 4;
174 if ((i == 0 and j == 1) or (i == 1 and j == 0)) // _01 and _10 --> 5
175 return 5;
176 assert(i < 3 and j < 3 && "For Voigt notation the indices need to be 0,1 or 2.");
177 __builtin_unreachable();
178 }
179
194 template <typename ScalarType = double>
195 Eigen::Matrix<ScalarType, 6, 6> toVoigt(const Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3, 3, 3>>& ft) {
196 Eigen::Matrix<ScalarType, 6, 6> mat;
197 for (Eigen::Index i = 0; i < 3; ++i)
198 for (Eigen::Index j = 0; j < 3; ++j)
199 for (Eigen::Index k = 0; k < 3; ++k)
200 for (Eigen::Index l = 0; l < 3; ++l)
201 mat(toVoigt(i, j), toVoigt(k, l)) = ft(i, j, k, l);
202 return mat;
203 }
204
223 template <typename ST, int size, int Options>
224 requires(size > 0 and size <= 3) auto toVoigt(const Eigen::Matrix<ST, size, size, Options, size, size>& E,
225 bool isStrain = true) {
226 Eigen::Vector<ST, (size * (size + 1)) / 2> EVoigt;
227 EVoigt.template segment<size>(0) = E.diagonal();
228
229 const ST possibleStrainFactor = isStrain ? 2.0 : 1.0;
230 if constexpr (size == 2)
231 EVoigt(2) = E(0, 1) * possibleStrainFactor;
232 else if constexpr (size == 3) {
233 EVoigt(size) = E(1, 2) * possibleStrainFactor;
234 EVoigt(size + 1) = E(0, 2) * possibleStrainFactor;
235 EVoigt(size + 2) = E(0, 1) * possibleStrainFactor;
236 }
237 return EVoigt;
238 }
239
255 template <typename ST, int size>
256 requires(size == 1 or size == 3 or size == 6) auto fromVoigt(const Eigen::Vector<ST, size>& EVoigt,
257 bool isStrain = true) {
258 constexpr int matrixSize = (-1 + ct_sqrt(1 + 8 * size)) / 2;
259 Eigen::Matrix<ST, matrixSize, matrixSize> E;
260 E.diagonal() = EVoigt.template segment<3>(0);
261
262 const ST possibleStrainFactor = isStrain ? 0.5 : 1.0;
263 if constexpr (matrixSize == 2)
264 E(0, 1) = E(1, 0) = EVoigt(2) * possibleStrainFactor;
265 else if constexpr (matrixSize == 3) {
266 E(2, 1) = E(1, 2) = EVoigt(matrixSize) * possibleStrainFactor;
267 E(2, 0) = E(0, 2) = EVoigt(matrixSize + 1) * possibleStrainFactor;
268 E(1, 0) = E(0, 1) = EVoigt(matrixSize + 2) * possibleStrainFactor;
269 }
270 return E;
271 }
272
286 constexpr std::array<size_t, 2> fromVoigt(size_t i) {
287 if (i < 3) // _00 -> 0, _11 -> 1, _22 -> 2
288 return {i, i};
289 else if (i == 3)
290 return {1, 2};
291 else if (i == 4)
292 return {0, 2};
293 else if (i == 5)
294 return {0, 1};
295 else {
296 assert(i < 6 && "For Voigt notation the indices need to be 0 and 5.");
297 __builtin_unreachable();
298 }
299 }
300
312 template <typename ScalarType>
313 auto fromVoigt(const Eigen::Matrix<ScalarType, 6, 6>& CVoigt) {
314 Eigen::TensorFixedSize<ScalarType, Eigen::Sizes<3, 3, 3, 3>> C;
315 // size_t iR=0,jR=0;
316 for (size_t i = 0; i < 6; ++i) {
317 for (size_t j = 0; j < 6; ++j) {
318 auto firstIndices = fromVoigt(i);
319 auto secondIndices = fromVoigt(j);
320 C(firstIndices[0], firstIndices[1], secondIndices[0], secondIndices[1]) = CVoigt(i, j);
321 }
322 }
323 return C;
324 }
325
326} // namespace Ikarus
Several concepts.
Implementation of math related algorithms.
auto fromVoigt(const Eigen::Vector< ST, size > &EVoigt, bool isStrain=true)
Converts a vector given in Voigt notation to a matrix.
Definition: tensorutils.hh:256
auto symmetricIdentityFourthOrder()
Generates a symmetric identity fourth-order tensor.
Definition: tensorutils.hh:60
auto symmetricFourthOrder(const auto &A, const auto &B)
Generates a symmetric fourth-order tensor based on two second-order input tensors.
Definition: tensorutils.hh:82
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:167
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:122
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:101
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:146
Definition: simpleassemblers.hh:21
constexpr T ct_sqrt(T x)
Compile-time square root for integer types.
Definition: math.hh:46
Definition: concepts.hh:25