version 0.4.1
resultevaluators.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
11#pragma once
12
13#include <dune/common/math.hh>
14
16
18
26{
34 template <typename R>
35 double operator()(const R& resultArray, [[maybe_unused]] const int comp) const {
36 const auto s_x = resultArray(0, 0);
37 const auto s_y = resultArray(1, 0);
38 if constexpr (R::CompileTimeTraits::RowsAtCompileTime == 3) {
39 const auto s_xy = resultArray(2, 0);
40 return std::sqrt(Dune::power(s_x, 2) + Dune::power(s_y, 2) - s_x * s_y + 3 * Dune::power(s_xy, 2));
41 } else {
42 const auto s_z = resultArray(2, 0);
43 const auto s_yz = resultArray(3, 0);
44 const auto s_xz = resultArray(4, 0);
45 const auto s_xy = resultArray(5, 0);
46
47 return std::sqrt(Dune::power(s_x, 2) + Dune::power(s_y, 2) + Dune::power(s_z, 2) - s_x * s_y - s_x * s_z -
48 s_y * s_z + 3 * (Dune::power(s_xy, 2) + Dune::power(s_xz, 2) + Dune::power(s_yz, 2)));
49 }
50 }
51
56 static std::string name() { return "VonMises"; }
57
62 static int ncomps() { return 1; }
63};
64
72template <int dim>
73requires(dim == 2 or dim == 3)
75{
82 double operator()(const auto& resultArray, const int comp) const {
83 auto mat = fromVoigt(resultArray, false);
84 Eigen::SelfAdjointEigenSolver<decltype(mat)> eigensolver(mat, Eigen::EigenvaluesOnly);
85 return eigensolver.eigenvalues()[dim - 1 - comp];
86 }
87
92 static std::string name() { return "PrincipalStress"; }
93
98 static int ncomps() { return dim; }
99};
100
101} // namespace Ikarus::ResultEvaluators
Helper for the Eigen::Tensor types.
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: resultevaluators.hh:17
Struct for calculating von Mises stress.
Definition: resultevaluators.hh:26
double operator()(const R &resultArray, const int comp) const
Calculate the result quantity (von Mises stress)
Definition: resultevaluators.hh:35
static std::string name()
Get the name of the result type (VonMises)
Definition: resultevaluators.hh:56
static int ncomps()
Get the number of components in the result (always 1 for VonMises)
Definition: resultevaluators.hh:62
Struct for calculating principal stresses.
Definition: resultevaluators.hh:75
double operator()(const auto &resultArray, const int comp) const
Calculate the result quantity (principal stress)
Definition: resultevaluators.hh:82
static std::string name()
Get the name of the result type (PrincipalStress)
Definition: resultevaluators.hh:92
static int ncomps()
Get the number of components in the result.
Definition: resultevaluators.hh:98