version 0.4.2
materialhelpers.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers ikarus@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#pragma once
5
12#include <ranges>
13
14#include <dune/common/float_cmp.hh>
15
16#include <Eigen/Dense>
17
20
21namespace Ikarus::Materials {
26{
27 Eigen::Index row;
28 Eigen::Index col;
29};
30} // namespace Ikarus::Materials
31namespace Ikarus::Impl {
32
39template <size_t size>
40consteval auto createfreeVoigtIndices(const std::array<Materials::MatrixIndexPair, size>& fixed) {
41 std::array<size_t, 6 - size> res{};
42 std::array<size_t, size> voigtFixedIndices;
43 std::ranges::transform(fixed, voigtFixedIndices.begin(), [](auto pair) { return toVoigt(pair.row, pair.col); });
44 std::ranges::sort(voigtFixedIndices);
45 std::ranges::set_difference(std::ranges::iota_view(size_t(0), size_t(6)), voigtFixedIndices, res.begin());
46 std::ranges::sort(res);
47 return res;
48}
49
56template <size_t size>
57consteval auto createFixedVoigtIndices(const std::array<Materials::MatrixIndexPair, size>& fixed) {
58 std::array<size_t, size> fixedIndices;
59 std::ranges::transform(fixed, fixedIndices.begin(), [](auto pair) { return toVoigt(pair.row, pair.col); });
60 std::ranges::sort(fixedIndices);
61 return fixedIndices;
62}
63
70template <size_t size>
71constexpr size_t countDiagonalIndices(const std::array<Materials::MatrixIndexPair, size>& fixed) {
72 size_t count = 0;
73 for (auto v : fixed) {
74 if (v.col == v.row)
75 ++count;
76 }
77 return count;
78}
79
88template <typename Derived>
89auto maybeFromVoigt(const Eigen::MatrixBase<Derived>& v, bool isStrain = true) {
90 if constexpr (Concepts::EigenVector<Derived>) { // receiving vector means Voigt notation
91 return fromVoigt(v.derived(), isStrain).eval();
92 } else
93 return v.derived().eval();
94}
95
104template <typename Derived>
105auto maybeToVoigt(const Eigen::MatrixBase<Derived>& v, bool isStrain = false) {
106 if constexpr (Concepts::EigenVector<Derived>) { // receiving vector means Voigt notation
107 return v.derived().eval();
108 } else
109 return toVoigt(v.derived(), isStrain).eval();
110}
111
119template <typename ScalarType>
120void checkPositiveOrAbort(ScalarType det, double eps = 1e-10) {
121 if (Dune::FloatCmp::le(static_cast<double>(det), 0.0, eps)) {
122 std::cerr << "Determinant of right Cauchy Green tensor C must be greater than zero. detC = " +
123 std::to_string(static_cast<double>(det));
124 abort();
125 }
126}
127
138template <typename ScalarType, typename Derived>
139auto principalStretches(const Eigen::MatrixBase<Derived>& C, int options = Eigen::ComputeEigenvectors) {
140 Eigen::SelfAdjointEigenSolver<Derived> eigensolver{};
141
142 eigensolver.compute(C, options);
143
144 if (eigensolver.info() != Eigen::Success)
145 DUNE_THROW(Dune::MathError, "Failed to compute eigenvalues and eigenvectors of C.");
146
147 auto& eigenvalues = eigensolver.eigenvalues();
148 auto& eigenvectors = options == Eigen::ComputeEigenvectors ? eigensolver.eigenvectors() : Derived::Zero();
149
150 auto principalStretches = eigenvalues.cwiseSqrt().eval();
151 return std::make_pair(principalStretches, eigenvectors);
152}
153
161template <Concepts::EigenVector3 Vector>
162inline Vector::Scalar determinantFromPrincipalValues(const Vector& principalValues) {
163 return principalValues.prod();
164}
165
174template <Concepts::EigenVector3 Vector>
175inline Vector deviatoricStretches(const Vector& lambda) {
176 using ScalarType = typename Vector::Scalar;
177 ScalarType J = determinantFromPrincipalValues(lambda);
178 ScalarType Jmod = pow(J, -1.0 / 3.0);
179 return Jmod * lambda;
180}
181
189template <Concepts::EigenVector3 Vector>
190inline Vector invariants(const Vector& lambda) {
191 using ScalarType = typename Vector::Scalar;
192 auto lambdaSquared = lambda.array().square();
193 auto invariants = Vector::Zero().eval();
194
195 invariants[0] = lambdaSquared.sum();
196 invariants[1] =
197 lambdaSquared[0] * lambdaSquared[1] + lambdaSquared[1] * lambdaSquared[2] + lambdaSquared[0] * lambdaSquared[2];
198 invariants[2] = determinantFromPrincipalValues(lambdaSquared);
199
200 return invariants;
201}
202} // namespace Ikarus::Impl
Helper for the Eigen::Tensor types.
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
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:182
Definition: arrudaboyce.hh:27
Represents a pair of stress or strain matrix indices (row and column).
Definition: materialhelpers.hh:26
Eigen::Index col
Column index.
Definition: materialhelpers.hh:28
Eigen::Index row
Row index.
Definition: materialhelpers.hh:27
Several concepts.