13#include <dune/common/fvector.hh>
14#include <dune/localfefunctions/eigenDuneTransformations.hh>
15#include <dune/localfefunctions/expressions/linearStrainsExpr.hh>
40 template <
typename GEO,
typename EAST>
42 const EAST& easFunction,
const auto& alpha) {
43 using namespace Dune::DerivativeDirections;
45 constexpr int myDim = GEO::mydimension;
46 const auto& strainFunction = Dune::linearStrains(uFunction);
47 const auto E = (strainFunction.evaluate(gpPos, on(gridElement))).eval();
48 const auto Etilde = (easFunction(gpPos) * alpha).eval();
49 return (E + Etilde).eval();
72 template <
int wrtCoeff,
typename GEO,
typename EAST>
73 static auto firstDerivative(
const GEO& geo,
const auto& uFunction,
const auto& localBasis,
const auto& gpIndex,
75 const auto& alpha,
const int node = sNaN) {
76 if constexpr (wrtCoeff == 0) {
77 using namespace Dune::DerivativeDirections;
79 constexpr int myDim = GEO::mydimension;
80 const auto& strainFunction = Dune::linearStrains(uFunction);
81 const auto bopI = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(node)), on(gridElement));
83 }
else if constexpr (wrtCoeff == 1) {
84 return easFunction(gpPos);
86 static_assert(Dune::AlwaysFalse<GEO>::value,
87 "firstDerivative can only be called with wrtCoeff as 0 and 1 indicating first derivative of the "
88 "linear strain w.r.t d and "
89 "alpha, respectively.");
116 template <
int wrtCoeff,
typename GEO,
typename ST,
typename EAST>
117 static auto secondDerivative(
const GEO& geo,
const auto& uFunction,
const auto& localBasis,
const auto& gpIndex,
119 const Eigen::Vector<ST, GEO::mydimension*(GEO::mydimension + 1) / 2>& S,
120 const EAST& easFunction,
const auto& alpha,
const int I = sNaN,
const int J = sNaN) {
121 constexpr int myDim = GEO::mydimension;
122 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
123 if constexpr (wrtCoeff == 0) {
124 const auto kg = Eigen::Matrix<ST, myDim, myDim>::Zero().eval();
126 }
else if constexpr (wrtCoeff == 1) {
127 const auto kg = Eigen::Matrix<ST, enhancedStrainSize, enhancedStrainSize>::Zero().eval();
129 }
else if constexpr (wrtCoeff == 2) {
130 const auto kg = Eigen::Matrix<ST, enhancedStrainSize, myDim>::Zero().eval();
134 Dune::AlwaysFalse<GEO>::value,
135 "secondDerivative can only be called with wrtCoeff as 0, 1 and 2 indicating second derivative of the "
136 "linear strain w.r.t d and "
137 "alpha and the mixed derivative, respectively.");
141 static constexpr auto name() {
return std::string(
"Linear Strain"); }
144 static constexpr int sNaN = std::numeric_limits<int>::signaling_NaN();
Definition: greenlagrangestrain.hh:19
Definition: utils/dirichletvalues.hh:30
A struct computing the value, first and second derivatives of the linear strain tensor (in Voigt nota...
Definition: linearstrain.hh:25
static auto value(const GEO &geo, const auto &uFunction, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const EAST &easFunction, const auto &alpha)
Compute the strain vector at a given integration point or its index.
Definition: linearstrain.hh:41
static auto firstDerivative(const GEO &geo, const auto &uFunction, const auto &localBasis, const auto &gpIndex, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const EAST &easFunction, const auto &alpha, const int node=sNaN)
Compute the first derivative of the linear strain w.r.t d or alpha for a given node and integration p...
Definition: linearstrain.hh:73
static constexpr auto name()
The name of the strain measure enhanced w.r.t EAS method.
Definition: linearstrain.hh:141
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 EAST &easFunction, const auto &alpha, const int I=sNaN, const int J=sNaN)
Compute the second derivative of the linear strain w.r.t d or alpha for a given node and integration ...
Definition: linearstrain.hh:117
Definition: utils/dirichletvalues.hh:32