version 0.4.4
pk2stress.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers ikarus@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
11#pragma once
12
13#include <dune/common/fvector.hh>
14#include <dune/localfefunctions/eigenDuneTransformations.hh>
15#include <dune/localfefunctions/expressions/greenLagrangeStrains.hh>
16
17#include <Eigen/Core>
18
20
21namespace Ikarus::PS {
22
24{
38 template <typename GEO, typename AST>
39 static auto value(const GEO& geo, const Dune::FieldVector<double, GEO::mydimension>& gpPos, const AST& asFunction,
40 const auto& beta) {
41 const auto P = asFunction(gpPos);
42 return (P * beta).eval();
43 }
44
64 template <typename GEO, typename AST>
65 static auto firstDerivative(const GEO& geo, const auto& uFunction, const auto& localBasis, const auto& gpIndex,
66 const Dune::FieldVector<double, GEO::mydimension>& gpPos, const AST& asFunction,
67 const auto& beta, const int node = sNaN) {
68 return asFunction(gpPos);
69 }
70
93 template <typename GEO, typename ST, typename AST>
94 static auto secondDerivative(const GEO& geo, const auto& uFunction, const auto& localBasis, const auto& gpIndex,
96 const Eigen::Vector<ST, GEO::mydimension*(GEO::mydimension + 1) / 2>& S,
97 const AST& asFunction, const auto& beta, const int I = sNaN, const int J = sNaN) {
98 constexpr int myDim = GEO::mydimension;
99 constexpr int assumedStressSize = AST::assumedStressSize;
100
101 return Eigen::Matrix<ST, assumedStressSize, assumedStressSize>::Zero().eval();
102 }
103
105 static constexpr auto name() { return std::string("PK2 Stress"); }
106
107private:
108 static constexpr int sNaN = std::numeric_limits<int>::signaling_NaN();
109};
110
111} // namespace Ikarus::PS
Helper for the Eigen::Tensor types.
Definition: linearstress.hh:19
Definition: pk2stress.hh:24
static constexpr auto name()
The name of the assumed stress type.
Definition: pk2stress.hh:105
static auto value(const GEO &geo, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const AST &asFunction, const auto &beta)
Compute the stress vector at a given integration point or its index.
Definition: pk2stress.hh:39
static auto secondDerivative(const GEO &geo, const auto &uFunction, const auto &localBasis, const auto &gpIndex, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const Eigen::Vector< ST, GEO::mydimension *(GEO::mydimension+1)/2 > &S, const AST &asFunction, const auto &beta, const int I=sNaN, const int J=sNaN)
Compute the second derivative of the PK2 stress w.r.t beta for a given node and integration point.
Definition: pk2stress.hh:94
static auto firstDerivative(const GEO &geo, const auto &uFunction, const auto &localBasis, const auto &gpIndex, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const AST &asFunction, const auto &beta, const int node=sNaN)
Compute the first derivative of the PK2 stress w.r.t beta for a given node and integration point.
Definition: pk2stress.hh:65
Definition: utils/dirichletvalues.hh:32