13#include <dune/common/math.hh>
19 template <
typename ScalarType,
int size>
36 template <
typename ElementType,
typename FERequirements,
int size,
typename ScalarType>
41 fe.calculateAt(req, pos, res_);
45 const auto s_x = sigma(0, 0);
46 const auto s_y = sigma(1, 0);
47 const auto s_xy = sigma(2, 0);
49 return std::sqrt(std::pow(s_x, 2) + Dune::power(s_y, 2) - s_x * s_y + 3 * Dune::power(s_xy, 2));
56 static std::string
name() {
return "VonMises"; }
77 template <
typename ElementType,
typename FERequirements,
int size,
typename ScalarType>
81 fe.calculateAt(req, pos, res_);
85 const auto s_x = sigma(0, 0);
86 const auto s_y = sigma(1, 0);
87 const auto s_xy = sigma(2, 0);
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));
92 return comp == 0 ? t1 + t2 : t1 - t2;
99 static std::string
name() {
return "PrincipalStress"; }
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