version 0.4
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
17namespace Dune {
18 // Forward declaration
19 template <typename ScalarType, int size>
21} // namespace Dune
22
24
35 struct VonMises {
36 template <typename ElementType, typename FERequirements, int size, typename ScalarType>
37 double operator()(const ElementType& fe, const ResultRequirements<FERequirements>& req,
38 const Dune::FieldVector<ScalarType, size>& pos, [[maybe_unused]] int comp) const
39 requires(size == 2) {
41 fe.calculateAt(req, pos, res_);
42
43 const auto& [resultType, sigma] = res_.getSingleResult();
44 assert(resultType == ResultType::cauchyStress or resultType == ResultType::PK2Stress);
45 const auto s_x = sigma(0, 0);
46 const auto s_y = sigma(1, 0);
47 const auto s_xy = sigma(2, 0);
48
49 return std::sqrt(std::pow(s_x, 2) + Dune::power(s_y, 2) - s_x * s_y + 3 * Dune::power(s_xy, 2));
50 }
51
56 static std::string name() { return "VonMises"; }
57
62 static int ncomps() { return 1; }
63 };
64
77 template <typename ElementType, typename FERequirements, int size, typename ScalarType>
78 double operator()(const ElementType& fe, const ResultRequirements<FERequirements>& req,
79 const Dune::FieldVector<ScalarType, size>& pos, int comp) const requires(size == 2) {
81 fe.calculateAt(req, pos, res_);
82
83 const auto& [resultType, sigma] = res_.getSingleResult();
84 assert(resultType == ResultType::cauchyStress or resultType == ResultType::PK2Stress);
85 const auto s_x = sigma(0, 0);
86 const auto s_y = sigma(1, 0);
87 const auto s_xy = sigma(2, 0);
88
89 auto t1 = (s_x + s_y) / 2;
90 auto t2 = std::sqrt(Dune::power((s_x - s_y) / 2, 2) + Dune::power(s_xy, 2));
91
92 return comp == 0 ? t1 + t2 : t1 - t2;
93 }
94
99 static std::string name() { return "PrincipalStress"; }
100
105 static int ncomps() { return 2; }
106 };
107
108} // namespace Ikarus::ResultEvaluators
Definition of the LinearElastic class for finite element mechanics computations.
Definition: resultevaluators.hh:17
Definition: resultevaluators.hh:23
Class representing a map of result types to result arrays.
Definition: ferequirements.hh:342
auto & getSingleResult()
Get the result array for a single result type.
Definition: ferequirements.hh:379
Class representing the requirements for obtaining specific results.
Definition: ferequirements.hh:405
Definition: resultevaluators.hh:20
Struct for calculating von Mises stress.
Definition: resultevaluators.hh:35
double operator()(const ElementType &fe, const ResultRequirements< FERequirements > &req, const Dune::FieldVector< ScalarType, size > &pos, int comp) const
Definition: resultevaluators.hh:37
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:76
double operator()(const ElementType &fe, const ResultRequirements< FERequirements > &req, const Dune::FieldVector< ScalarType, size > &pos, int comp) const
Definition: resultevaluators.hh:78
static std::string name()
Get the name of the result type (PrincipalStress)
Definition: resultevaluators.hh:99
static int ncomps()
Get the number of components in the result (always 2 for PrincipalStress)
Definition: resultevaluators.hh:105