version 0.4.1
greenlagrangestrain.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 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/fvector.hh>
14#include <dune/localfefunctions/eigenDuneTransformations.hh>
15#include <dune/localfefunctions/expressions/greenLagrangeStrains.hh>
16
17#include <Eigen/Core>
18
19namespace Ikarus::EAS {
25{
40 template <typename GEO, typename EAST>
41 static auto value(const GEO& geo, const auto& uFunction, const Dune::FieldVector<double, GEO::mydimension>& gpPos,
42 const EAST& easFunction, const auto& alpha) {
43 using namespace Dune::DerivativeDirections;
44 using namespace Dune;
45 constexpr int myDim = GEO::mydimension;
46 const auto& strainFunction = Dune::greenLagrangeStrains(uFunction);
47 const auto E = (strainFunction.evaluate(gpPos, on(gridElement))).eval();
48 const auto Etilde = (easFunction(gpPos) * alpha).eval();
49 return (E + Etilde).eval();
50 }
51
73 template <int wrtCoeff, typename GEO, typename EAST>
74 static auto firstDerivative(const GEO& geo, const auto& uFunction, const auto& localBasis, const auto& gpIndex,
75 const Dune::FieldVector<double, GEO::mydimension>& gpPos, const EAST& easFunction,
76 const auto& alpha, const int node = sNaN) {
77 if constexpr (wrtCoeff == 0) {
78 using namespace Dune::DerivativeDirections;
79 using namespace Dune;
80 constexpr int myDim = GEO::mydimension;
81 const auto& strainFunction = Dune::greenLagrangeStrains(uFunction);
82 const auto bopI = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(node)), on(gridElement));
83 return bopI.eval();
84 } else if constexpr (wrtCoeff == 1) {
85 return easFunction(gpPos);
86 } else
87 static_assert(Dune::AlwaysFalse<GEO>::value,
88 "firstDerivative can only be called with wrtCoeff as 0 and 1 indicating first derivative of the "
89 "Green-Lagrange strain w.r.t d and "
90 "alpha, respectively.");
91 }
92
118 template <int wrtCoeff, typename GEO, typename ST, typename EAST>
119 static auto secondDerivative(const GEO& geo, const auto& uFunction, const auto& localBasis, const auto& gpIndex,
121 const Eigen::Vector<ST, GEO::mydimension*(GEO::mydimension + 1) / 2>& S,
122 const EAST& easFunction, const auto& alpha, const int I = sNaN, const int J = sNaN) {
123 constexpr int myDim = GEO::mydimension;
124 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
125 if constexpr (wrtCoeff == 0) {
126 using namespace Dune::DerivativeDirections;
127 using namespace Dune;
128 const auto& strainFunction = Dune::greenLagrangeStrains(uFunction);
129 const auto kg = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(I, J)), along(S), on(gridElement));
130 return kg.eval(); // E_dd
131 } else if constexpr (wrtCoeff == 1) {
132 const auto kg = Eigen::Matrix<ST, enhancedStrainSize, enhancedStrainSize>::Zero().eval();
133 return kg; // E_aa
134 } else if constexpr (wrtCoeff == 2) {
135 const auto kg = Eigen::Matrix<ST, enhancedStrainSize, myDim>::Zero().eval();
136 return kg; // E_ad
137 } else
138 static_assert(
139 Dune::AlwaysFalse<GEO>::value,
140 "secondDerivative can only be called with wrtCoeff as 0, 1 and 2 indicating second derivative of the "
141 "Green-Lagrange strain w.r.t d and "
142 "alpha and the mixed derivative, respectively.");
143 }
144
146 static constexpr auto name() { return std::string("Green-Lagrange Strain"); }
147
148private:
149 static constexpr int sNaN = std::numeric_limits<int>::signaling_NaN();
150};
151
152} // namespace Ikarus::EAS
Definition: easfunctions/displacementgradient.hh:21
Definition: utils/dirichletvalues.hh:30
A struct computing the value, first and second derivatives of the Green-Lagrange strain tensor (in Vo...
Definition: greenlagrangestrain.hh:25
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 Green-Lagrange strain w.r.t d or alpha for a given node and inte...
Definition: greenlagrangestrain.hh:119
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: greenlagrangestrain.hh:41
static constexpr auto name()
The name of the strain measure enhanced w.r.t EAS method.
Definition: greenlagrangestrain.hh:146
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 Green-Lagrange strain w.r.t d or alpha for a given node and integ...
Definition: greenlagrangestrain.hh:74
Definition: utils/dirichletvalues.hh:32