version 0.4.2
vanishingstress.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
10// SPDX-License-Identifier: LGPL-3.0-or-later
11
12#pragma once
13
18
19namespace Ikarus {
20
21namespace Impl {
22
26 struct StressIndexPair
27 {
28 Eigen::Index row;
29 Eigen::Index col;
30 };
31
38 template <size_t size>
39 consteval auto createfreeVoigtIndices(const std::array<StressIndexPair, size>& fixed) {
40 std::array<size_t, 6 - size> res{};
41 std::array<size_t, size> voigtFixedIndices;
42 std::ranges::transform(fixed, voigtFixedIndices.begin(), [](auto pair) { return toVoigt(pair.row, pair.col); });
43 std::ranges::sort(voigtFixedIndices);
44 std::ranges::set_difference(std::ranges::iota_view(size_t(0), size_t(6)), voigtFixedIndices, res.begin());
45 std::ranges::sort(res);
46 return res;
47 }
48
55 template <size_t size>
56 consteval auto createFixedVoigtIndices(const std::array<StressIndexPair, size>& fixed) {
57 std::array<size_t, size> fixedIndices;
58 std::ranges::transform(fixed, fixedIndices.begin(), [](auto pair) { return toVoigt(pair.row, pair.col); });
59 std::ranges::sort(fixedIndices);
60 return fixedIndices;
61 }
62
69 template <size_t size>
70 constexpr size_t countDiagonalIndices(const std::array<StressIndexPair, size>& fixed) {
71 size_t count = 0;
72 for (auto v : fixed) {
73 if (v.col == v.row)
74 ++count;
75 }
76 return count;
77 }
78
79} // namespace Impl
80
87template <auto stressIndexPair, typename MI>
88struct VanishingStress : public Material<VanishingStress<stressIndexPair, MI>>
89{
95 explicit VanishingStress(MI mat, typename MI::ScalarType tol = 1e-12)
96 : matImpl_{mat},
97 tol_{tol} {}
98
99 using Underlying = MI;
100
101 static constexpr auto fixedPairs = stressIndexPair;
102 static constexpr auto freeVoigtIndices = createfreeVoigtIndices(fixedPairs);
103 static constexpr auto fixedVoigtIndices = createFixedVoigtIndices(fixedPairs);
104 static constexpr auto fixedDiagonalVoigtIndicesSize =
105 countDiagonalIndices(fixedPairs);
106 static constexpr auto freeStrains = freeVoigtIndices.size();
107 using ScalarType = typename Underlying::ScalarType;
108
109 [[nodiscard]] constexpr std::string nameImpl() const noexcept {
110 auto matName = matImpl_.name() + "_Vanishing(";
111 for (auto p : fixedPairs)
112 matName += "(" + std::to_string(p.row) + std::to_string(p.col) + ")";
113 matName += ")";
114 return matName;
115 }
116
117 static constexpr auto strainTag = Underlying::strainTag;
118 static constexpr auto stressTag = Underlying::stressTag;
119 static constexpr auto tangentModuliTag = Underlying::tangentModuliTag;
120 static constexpr bool energyAcceptsVoigt = Underlying::energyAcceptsVoigt;
121 static constexpr bool stressToVoigt = true;
122 static constexpr bool stressAcceptsVoigt = true;
123 static constexpr bool moduliToVoigt = true;
124 static constexpr bool moduliAcceptsVoigt = true;
125 static constexpr double derivativeFactor = 1;
126
133 template <typename Derived>
134 ScalarType storedEnergyImpl(const Eigen::MatrixBase<Derived>& E) const {
135 const auto [nonOp, Esol] = reduceStress(E);
136 return matImpl_.storedEnergyImpl(Esol);
137 }
138
146 template <bool voigt, typename Derived>
147 auto stressesImpl(const Eigen::MatrixBase<Derived>& E) const {
148 const auto [nonOp, Esol] = reduceStress(E);
149 auto stressesRed = matImpl_.template stresses<Underlying::strainTag, true>(Esol);
150
151 if constexpr (voigt) {
152 return removeCol(stressesRed, fixedVoigtIndices);
153 } else
154 return fromVoigt(stressesRed, false);
155 }
156
164 template <bool voigt, typename Derived>
165 auto tangentModuliImpl(const Eigen::MatrixBase<Derived>& E) const {
166 const auto [nonOp, Esol] = reduceStress(E);
167 auto C = matImpl_.template tangentModuli<Underlying::strainTag, true>(Esol);
168 if constexpr (voigt)
170 else
171 return fromVoigt(C);
172 }
173
179 template <typename ScalarTypeOther>
180 auto rebind() const {
181 auto reboundMatImpl = matImpl_.template rebind<ScalarTypeOther>();
183 }
184
185private:
192 template <typename Derived>
193 decltype(auto) maybeFromVoigt(const Eigen::MatrixBase<Derived>& E) const {
194 if constexpr (Concepts::EigenVector<Derived>) { // receiving vector means Voigt notation
195 return fromVoigt(E.derived(), true);
196 } else
197 return E.derived();
198 }
199
205 template <typename Derived>
206 void initUnknownStrains(Eigen::MatrixBase<Derived>& E) const {
207 for (size_t i = 0; i < fixedPairs.size(); ++i) {
208 ScalarType initialVal = E(fixedPairs[i].row, fixedPairs[i].col);
210 if (Dune::FloatCmp::eq(initialVal, ScalarType(0.0)) and (fixedPairs[i].row == fixedPairs[i].col))
211 initialVal = ScalarType(1.0);
212 }
213 if (fixedPairs[i].row != fixedPairs[i].col)
214 initialVal = ScalarType(0.0);
215 E(fixedPairs[i].row, fixedPairs[i].col) = E(fixedPairs[i].col, fixedPairs[i].row) = initialVal;
216 }
217 }
218
225 template <typename Derived>
226 auto reduceStress(const Eigen::MatrixBase<Derived>& Eraw) const {
227 auto E = maybeFromVoigt(Eraw);
228 initUnknownStrains(E);
229
230 std::array<size_t, fixedDiagonalVoigtIndicesSize> fixedDiagonalVoigtIndices;
231 for (size_t ri = 0; auto i : fixedVoigtIndices) {
232 auto indexPair = fromVoigt(i);
233 if (indexPair[0] == indexPair[1])
234 fixedDiagonalVoigtIndices[ri++] = i;
235 }
236
237 auto f = [&](auto&) {
238 auto S = matImpl_.template stresses<Underlying::strainTag, true>(E);
239 return S(fixedDiagonalVoigtIndices).eval();
240 };
241 auto df = [&](auto&) {
242 auto moduli = (matImpl_.template tangentModuli<Underlying::strainTag, true>(E)).eval();
243 return (moduli(fixedDiagonalVoigtIndices, fixedDiagonalVoigtIndices) / Underlying::derivativeFactor).eval();
244 };
245
246 auto Er = E(fixedDiagonalVoigtIndices, fixedDiagonalVoigtIndices).eval().template cast<ScalarType>();
247 auto nonOp = Ikarus::NonLinearOperator(functions(f, df), parameter(Er));
248
249 auto linearSolver = [](auto& r, auto& A) { return (A.inverse() * r).eval(); };
250 auto updateFunction = [&](auto& /* Ex33 */, auto& ecomps) {
251 for (int ri = 0; auto i : fixedDiagonalVoigtIndices) {
252 auto indexPair = fromVoigt(i);
253 E(indexPair[0], indexPair[1]) += ecomps(ri++);
254 }
255 };
256 // THE CTAD is broken for designated initializers in clang 16, when we drop support this can be simplified
257 NewtonRaphsonConfig<decltype(linearSolver), decltype(updateFunction)> nrs{
258 .parameters = {.tol = tol_, .maxIter = 100},
259 .linearSolver = linearSolver, .updateFunction = updateFunction
260 };
261
262 auto nr = createNonlinearSolver(std::move(nrs), nonOp);
263 if (!static_cast<bool>(nr->solve()))
264 DUNE_THROW(Dune::MathError, "The stress reduction of material " << nameImpl() << " was unsuccessful\n"
265 << "The strains are\n"
266 << E << "\n The stresses are\n"
267 << f(Er));
268 return std::make_pair(nonOp, E);
269 }
270
271 Underlying matImpl_;
272 double tol_{};
273};
274
283template <Impl::StressIndexPair... stressIndexPair, typename MaterialImpl>
284auto makeVanishingStress(MaterialImpl mat, typename MaterialImpl::ScalarType p_tol = 1e-12) {
285 return VanishingStress<std::to_array({stressIndexPair...}), MaterialImpl>(mat, p_tol);
286}
287
295template <typename MaterialImpl>
296auto planeStress(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
297 return makeVanishingStress<Impl::StressIndexPair{2, 1}, Impl::StressIndexPair{2, 0}, Impl::StressIndexPair{2, 2}>(
298 mat, tol);
299}
300
309template <typename MaterialImpl>
310auto shellMaterial(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
311 return makeVanishingStress<Impl::StressIndexPair{2, 2}>(mat, tol);
312}
313
322template <typename MaterialImpl>
323auto beamMaterial(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
324 return makeVanishingStress<Impl::StressIndexPair{1, 1}, Impl::StressIndexPair{2, 2}>(mat, tol);
325}
326} // namespace Ikarus
Provides a NonLinearOperator class for handling nonlinear operators.
Contains the generic NonlinearSolverFactory class.
Implementation of the Newton-Raphson method for solving nonlinear equations.
Contains the Material interface class and related template functions for material properties.
auto staticCondensation(const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfCondensedIndices > &indices)
Performs static condensation on a square matrix.
Definition: linearalgebrahelper.hh:498
auto removeCol(const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfRemovedCols > &indices)
Removes specified columns from a matrix.
Definition: linearalgebrahelper.hh:526
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
Definition: dirichletbcenforcement.hh:6
auto makeVanishingStress(MaterialImpl mat, typename MaterialImpl::ScalarType p_tol=1e-12)
Factory function to create a VanishingStress material with specified stress indices.
Definition: vanishingstress.hh:284
::value auto createNonlinearSolver(NRConfig &&config, NLO &&nonLinearOperator)
Function to create a NewtonRaphson solver instance.
Definition: newtonraphson.hh:64
auto shellMaterial(const MaterialImpl &mat, typename MaterialImpl::ScalarType tol=1e-8)
Factory function to create a VanishingStress material for a shell material with zero normal stress co...
Definition: vanishingstress.hh:310
auto beamMaterial(const MaterialImpl &mat, typename MaterialImpl::ScalarType tol=1e-8)
Factory function to create a VanishingStress material for a beam material with two zero normal stress...
Definition: vanishingstress.hh:323
auto functions(Args &&... args)
Creates a Functions object.
Definition: nonlinearoperator.hh:127
NonLinearOperator(const Impl::Functions< DerivativeArgs &&... > &a, const Impl::Parameter< ParameterArgs... > &b) -> NonLinearOperator< Impl::Functions< DerivativeArgs... >, Impl::Parameter< ParameterArgs... > >
auto parameter(Args &&... args)
Creates a Parameter object.
Definition: nonlinearoperator.hh:115
auto planeStress(const MaterialImpl &mat, typename MaterialImpl::ScalarType tol=1e-8)
Factory function to create a VanishingStress material for plane stress conditions.
Definition: vanishingstress.hh:296
Interface classf or materials.
Definition: interface.hh:77
VanishingStress material model that enforces stress components to be zero.
Definition: vanishingstress.hh:89
static constexpr auto tangentModuliTag
Tangent moduli tag.
Definition: vanishingstress.hh:119
auto rebind() const
Rebinds the material to a different scalar type.
Definition: vanishingstress.hh:180
static constexpr auto fixedPairs
Array of fixed stress components.
Definition: vanishingstress.hh:101
static constexpr double derivativeFactor
Derivative factor.
Definition: vanishingstress.hh:125
VanishingStress(MI mat, typename MI::ScalarType tol=1e-12)
Constructor for VanishingStress.
Definition: vanishingstress.hh:95
static constexpr auto strainTag
Strain tag.
Definition: vanishingstress.hh:117
static constexpr auto freeVoigtIndices
Free Voigt indices.
Definition: vanishingstress.hh:102
constexpr std::string nameImpl() const noexcept
Definition: vanishingstress.hh:109
static constexpr bool energyAcceptsVoigt
Energy accepts Voigt notation.
Definition: vanishingstress.hh:120
ScalarType storedEnergyImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the stored energy for the VanishingStress material.
Definition: vanishingstress.hh:134
MI Underlying
The underlying material type.
Definition: vanishingstress.hh:99
auto stressesImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the stresses for the VanishingStress material.
Definition: vanishingstress.hh:147
typename Underlying::ScalarType ScalarType
Scalar type.
Definition: vanishingstress.hh:107
static constexpr bool moduliToVoigt
Moduli to Voigt notation.
Definition: vanishingstress.hh:123
static constexpr bool stressToVoigt
Stress to Voigt notation.
Definition: vanishingstress.hh:121
auto tangentModuliImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the tangent moduli for the VanishingStress material.
Definition: vanishingstress.hh:165
static constexpr auto fixedDiagonalVoigtIndicesSize
Number of fixed diagonal indices.
Definition: vanishingstress.hh:104
static constexpr auto freeStrains
Number of free strains.
Definition: vanishingstress.hh:106
static constexpr auto stressTag
Stress tag.
Definition: vanishingstress.hh:118
static constexpr bool stressAcceptsVoigt
Stress accepts Voigt notation.
Definition: vanishingstress.hh:122
static constexpr auto fixedVoigtIndices
Fixed Voigt indices.
Definition: vanishingstress.hh:103
static constexpr bool moduliAcceptsVoigt
Moduli accepts Voigt notation.
Definition: vanishingstress.hh:124
Concept defining the requirements for Eigen vectors.
Definition: concepts.hh:352