version 0.4.1
vanishingstress.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
10#pragma once
11
12#include "materialhelpers.hh"
13
20
21namespace Ikarus::Materials {
22
29template <auto stressIndexPair, typename MI>
30struct VanishingStress : public Material<VanishingStress<stressIndexPair, MI>>
31{
32 using Underlying = MI;
33 using MaterialParameters = typename Underlying::MaterialParameters;
34 using StrainMatrix = typename Underlying::StrainMatrix;
35 using StressMatrix = typename Underlying::StressMatrix;
36 using MaterialTensor = typename Underlying::MaterialTensor;
37 static constexpr int dim = Underlying::dim;
38
39 static constexpr auto fixedPairs = stressIndexPair;
40 static constexpr auto freeVoigtIndices = Impl::createfreeVoigtIndices(fixedPairs);
41 static constexpr auto fixedVoigtIndices = Impl::createFixedVoigtIndices(fixedPairs);
42 static constexpr auto fixedDiagonalVoigtIndicesSize =
43 Impl::countDiagonalIndices(fixedPairs);
44 static constexpr auto freeStrains = freeVoigtIndices.size();
45 using ScalarType = typename Underlying::ScalarType;
47
48 static constexpr auto strainTag = Underlying::strainTag;
49 static constexpr auto stressTag = Underlying::stressTag;
50 static constexpr auto tangentModuliTag = Underlying::tangentModuliTag;
51 static constexpr bool energyAcceptsVoigt = Underlying::energyAcceptsVoigt;
52 static constexpr bool stressToVoigt = true;
53 static constexpr bool stressAcceptsVoigt = true;
54 static constexpr bool moduliToVoigt = true;
55 static constexpr bool moduliAcceptsVoigt = true;
56 static constexpr double derivativeFactorImpl = Underlying::derivativeFactorImpl;
57
63 explicit VanishingStress(MI mat, typename MI::ScalarType tol = 1e-12)
64 : matImpl_{mat},
65 tol_{tol} {}
66
67 [[nodiscard]] constexpr static std::string nameImpl() noexcept {
68 auto matName = MI::name() + "_VanishingStress(";
69 for (auto p : fixedPairs)
70 matName += "(" + std::to_string(p.row) + std::to_string(p.col) + ")";
71 matName += ")";
72 return matName;
73 }
74
78 MaterialParameters materialParametersImpl() const { return matImpl_.materialParametersImpl(); }
79
86 template <typename Derived>
87 ScalarType storedEnergyImpl(const Eigen::MatrixBase<Derived>& E) const {
88 const auto [nonOp, Esol] = reduceStress(E);
89 return matImpl_.storedEnergyImpl(Esol);
90 }
91
99 template <bool voigt, typename Derived>
100 auto stressesImpl(const Eigen::MatrixBase<Derived>& E) const {
101 const auto [nonOp, Esol] = reduceStress(E);
102 auto stressesRed = matImpl_.template stresses<Underlying::strainTag, true>(Esol);
103
104 if constexpr (voigt) {
105 return removeCol(stressesRed, fixedVoigtIndices);
106 } else
107 return fromVoigt(stressesRed, false);
108 }
109
117 template <bool voigt, typename Derived>
118 auto tangentModuliImpl(const Eigen::MatrixBase<Derived>& E) const {
119 const auto [nonOp, Esol] = reduceStress(E);
120 auto C = matImpl_.template tangentModuli<Underlying::strainTag, true>(Esol);
121 if constexpr (voigt)
123 else
124 return fromVoigt(C);
125 }
126
132 template <typename ScalarTypeOther>
133 auto rebind() const {
134 auto reboundMatImpl = matImpl_.template rebind<ScalarTypeOther>();
136 }
137
141 auto& underlying() const { return matImpl_; }
142
143private:
149 template <typename Derived>
150 void initUnknownStrains(Eigen::MatrixBase<Derived>& E) const {
151 for (size_t i = 0; i < fixedPairs.size(); ++i) {
152 ScalarType initialVal = E(fixedPairs[i].row, fixedPairs[i].col);
154 if (Dune::FloatCmp::eq(initialVal, ScalarType(0.0)) and (fixedPairs[i].row == fixedPairs[i].col))
155 initialVal = ScalarType(1.0);
156 }
157 if (fixedPairs[i].row != fixedPairs[i].col)
158 initialVal = ScalarType(0.0);
159 E(fixedPairs[i].row, fixedPairs[i].col) = E(fixedPairs[i].col, fixedPairs[i].row) = initialVal;
160 }
161 }
162
169 template <typename Derived>
170 auto reduceStress(const Eigen::MatrixBase<Derived>& Eraw) const {
171 auto E = Impl::maybeFromVoigt(Eraw);
172 initUnknownStrains(E);
173
174 std::array<size_t, fixedDiagonalVoigtIndicesSize> fixedDiagonalVoigtIndices;
175 for (size_t ri = 0; auto i : fixedVoigtIndices) {
176 auto indexPair = fromVoigt(i);
177 if (indexPair[0] == indexPair[1])
178 fixedDiagonalVoigtIndices[ri++] = i;
179 }
180
181 auto f = [&](const auto&) {
182 auto S = matImpl_.template stresses<Underlying::strainTag, true>(E);
183 return S(fixedDiagonalVoigtIndices).eval();
184 };
185 auto df = [&](const auto&) {
186 auto moduli = (matImpl_.template tangentModuli<Underlying::strainTag, true>(E)).eval();
187 return (moduli(fixedDiagonalVoigtIndices, fixedDiagonalVoigtIndices) / Underlying::derivativeFactor).eval();
188 };
189
190 auto Er = E(fixedDiagonalVoigtIndices, fixedDiagonalVoigtIndices).eval().template cast<ScalarType>();
191
192 Er.setZero();
193 auto diffFunction = Ikarus::makeDifferentiableFunction(functions(f, df), Er);
194
195 auto linearSolver = [](auto& r, auto& A) { return (A.inverse() * r).eval(); };
196 auto updateFunction = [&](auto&, const auto& ecomps) {
197 for (int ri = 0; auto i : fixedDiagonalVoigtIndices) {
198 auto indexPair = fromVoigt(i);
199 E(indexPair[0], indexPair[1]) += ecomps(ri++);
200 }
201 };
202
203 int minIter = isAutoDiff ? 1 : 0;
204 NewtonRaphsonConfig nrs({.tol = tol_, .maxIter = 100, .minIter = minIter}, linearSolver, updateFunction);
205
206 auto nr = createNonlinearSolver(std::move(nrs), diffFunction);
207 if (!static_cast<bool>(nr->solve(Er)))
208 DUNE_THROW(Dune::MathError, "The stress reduction of material " << nameImpl() << " was unsuccessful\n"
209 << "The strains are\n"
210 << E << "\n The stresses are\n"
211 << f(Er));
212 return std::make_pair(diffFunction, E);
213 }
214
215 Underlying matImpl_;
216 double tol_{};
217};
218
227template <MatrixIndexPair... stressIndexPair, typename MaterialImpl>
228auto makeVanishingStress(MaterialImpl mat, typename MaterialImpl::ScalarType p_tol = 1e-12) {
229 return VanishingStress<std::to_array({stressIndexPair...}), MaterialImpl>(mat, p_tol);
230}
231
239template <typename MaterialImpl>
240auto planeStress(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
241 return makeVanishingStress<MatrixIndexPair{2, 1}, MatrixIndexPair{2, 0}, MatrixIndexPair{2, 2}>(mat, tol);
242}
243
252template <typename MaterialImpl>
253auto shellMaterial(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
254 return makeVanishingStress<MatrixIndexPair{2, 2}>(mat, tol);
255}
256
265template <typename MaterialImpl>
266auto beamMaterial(const MaterialImpl& mat, typename MaterialImpl::ScalarType tol = 1e-8) {
267 return makeVanishingStress<MatrixIndexPair{1, 1}, MatrixIndexPair{2, 2}>(mat, tol);
268}
269} // namespace Ikarus::Materials
Helper for the Eigen::Tensor types.
Provides a DifferentiableFunction class for handling differentiable Functions.
Contains the generic NonlinearSolverFactory class.
Implementation of the Newton-Raphson method for solving nonlinear equations.
helper functions used by material model implementations.
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:539
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
::value auto createNonlinearSolver(NRConfig &&config, F &&f)
Function to create a NewtonRaphson solver instance.
Definition: newtonraphson.hh:82
auto makeDifferentiableFunction(const Impl::Functions< F... > &derivativesFunctions, const Arg &parameter)
Factory method for DifferentiableFunction It is a function taking several callables and the argument ...
Definition: differentiablefunction.hh:99
auto functions(Args &&... args)
Creates a Functions object.
Definition: differentiablefunction.hh:44
Definition: arrudaboyce.hh:27
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:228
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:253
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:240
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:266
Interface classf or materials.
Definition: finiteelements/mechanics/materials/interface.hh:80
Represents a pair of stress or strain matrix indices (row and column).
Definition: materialhelpers.hh:26
VanishingStress material model that enforces stress components to be zero.
Definition: vanishingstress.hh:31
auto tangentModuliImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the tangent moduli for the VanishingStress material.
Definition: vanishingstress.hh:118
static constexpr auto freeStrains
Number of free strains.
Definition: vanishingstress.hh:44
static constexpr auto freeVoigtIndices
Free Voigt indices.
Definition: vanishingstress.hh:40
static constexpr int dim
Definition: vanishingstress.hh:37
static constexpr bool moduliToVoigt
Moduli to Voigt notation.
Definition: vanishingstress.hh:54
typename Underlying::ScalarType ScalarType
Scalar type.
Definition: vanishingstress.hh:45
static constexpr auto fixedVoigtIndices
Fixed Voigt indices.
Definition: vanishingstress.hh:41
static constexpr bool moduliAcceptsVoigt
Moduli accepts Voigt notation.
Definition: vanishingstress.hh:55
typename Underlying::MaterialTensor MaterialTensor
Definition: vanishingstress.hh:36
static constexpr auto strainTag
Strain tag.
Definition: vanishingstress.hh:48
static constexpr auto stressTag
Stress tag.
Definition: vanishingstress.hh:49
MI Underlying
The underlying material type.
Definition: vanishingstress.hh:32
static constexpr bool isAutoDiff
Definition: vanishingstress.hh:46
auto stressesImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the stresses for the VanishingStress material.
Definition: vanishingstress.hh:100
static constexpr std::string nameImpl() noexcept
Definition: vanishingstress.hh:67
typename Underlying::MaterialParameters MaterialParameters
Definition: vanishingstress.hh:33
auto rebind() const
Rebinds the material to a different scalar type.
Definition: vanishingstress.hh:133
static constexpr double derivativeFactorImpl
Derivative factor.
Definition: vanishingstress.hh:56
ScalarType storedEnergyImpl(const Eigen::MatrixBase< Derived > &E) const
Computes the stored energy for the VanishingStress material.
Definition: vanishingstress.hh:87
auto & underlying() const
Returns a const reference to the underlying material.
Definition: vanishingstress.hh:141
static constexpr bool energyAcceptsVoigt
Energy accepts Voigt notation.
Definition: vanishingstress.hh:51
typename Underlying::StressMatrix StressMatrix
Definition: vanishingstress.hh:35
VanishingStress(MI mat, typename MI::ScalarType tol=1e-12)
Constructor for VanishingStress.
Definition: vanishingstress.hh:63
MaterialParameters materialParametersImpl() const
Returns the material parameters stored in the material.
Definition: vanishingstress.hh:78
static constexpr auto fixedPairs
Array of fixed stress components.
Definition: vanishingstress.hh:39
static constexpr bool stressToVoigt
Stress to Voigt notation.
Definition: vanishingstress.hh:52
static constexpr auto fixedDiagonalVoigtIndicesSize
Number of fixed diagonal indices.
Definition: vanishingstress.hh:42
static constexpr bool stressAcceptsVoigt
Stress accepts Voigt notation.
Definition: vanishingstress.hh:53
typename Underlying::StrainMatrix StrainMatrix
Definition: vanishingstress.hh:34
static constexpr auto tangentModuliTag
Tangent moduli tag.
Definition: vanishingstress.hh:50
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:615
Contains the Material interface class and related template functions for material properties.
Several concepts.