13#include <dune/common/fvector.hh>
14#include <dune/localfefunctions/eigenDuneTransformations.hh>
15#include <dune/localfefunctions/expressions/greenLagrangeStrains.hh>
44 template <
typename GEO,
typename EAST>
46 const EAST& easFunction,
const auto& alpha) {
47 using ST =
typename std::remove_cvref_t<
decltype(uFunction)>::ctype;
48 constexpr int myDim = GEO::mydimension;
49 const Eigen::Matrix<ST, myDim, myDim> H = computeDisplacementGradient(geo, uFunction, gpPos, easFunction, alpha);
50 const Eigen::Matrix<ST, myDim, myDim> E = 0.5 * (H + H.transpose() + H.transpose() * H);
74 template <
int wrtCoeff,
typename GEO,
typename EAST>
75 static auto firstDerivative(
const GEO& geo,
const auto& uFunction,
const auto& localBasis,
const auto& gpIndex,
77 const auto& alpha,
const int node = sNaN) {
78 using namespace Dune::DerivativeDirections;
80 using ST =
typename std::remove_cvref_t<
decltype(uFunction)>::ctype;
81 constexpr int myDim = GEO::mydimension;
82 static_assert(myDim == 2 or myDim == 3,
83 "An enhancement of displacement gradient can only be performed for the 2D or 3D case.");
84 constexpr int strainSize = myDim * (myDim + 1) / 2;
85 using MatrixType = Eigen::Matrix<ST, myDim, myDim>;
87 const MatrixType F = computeDeformationGradient(geo, uFunction, gpPos, easFunction, alpha)
90 const auto g0 = F.row(0);
91 const auto g1 = F.row(1);
93 if constexpr (wrtCoeff == 0) {
94 const auto gradUdI = uFunction.evaluateDerivative(gpIndex, wrt(spatialAll, coeff(node)), on(gridElement));
96 Eigen::Vector<ST, myDim> dN;
97 for (
const auto i : Dune::range(myDim))
98 dN[i] = getDiagonalEntry(gradUdI[i], 0);
99 Eigen::Matrix<ST, strainSize, myDim> bopI;
100 bopI.row(0) = dN[0] * g0;
101 bopI.row(1) = dN[1] * g1;
102 if constexpr (myDim == 2)
103 bopI.row(2) = dN[1] * g0 + dN[0] * g1;
105 const auto g2 = F.row(2);
106 bopI.row(2) = dN[2] * g2;
107 bopI.row(3) = dN[2] * g1 + dN[1] * g2;
108 bopI.row(4) = dN[2] * g0 + dN[0] * g2;
109 bopI.row(5) = dN[1] * g0 + dN[0] * g1;
112 }
else if constexpr (wrtCoeff == 1) {
113 const typename EAST::HType Harray = easFunction(gpPos);
114 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
115 Eigen::Matrix<ST, strainSize, enhancedStrainSize> M;
117 for (
const auto p : Dune::range(enhancedStrainSize)) {
118 const MatrixType Htilde = Harray[p];
119 M(0, p) = Htilde.col(0).dot(g0);
120 M(1, p) = Htilde.col(1).dot(g1);
121 if constexpr (myDim == 2)
122 M(2, p) = Htilde.col(1).dot(g0) + Htilde.col(0).dot(g1);
124 const auto g2 = F.row(2);
125 M(2, p) = Htilde.col(2).dot(g2);
126 M(3, p) = Htilde.col(2).dot(g1) + Htilde.col(1).dot(g2);
127 M(4, p) = Htilde.col(2).dot(g0) + Htilde.col(0).dot(g2);
128 M(5, p) = Htilde.col(1).dot(g0) + Htilde.col(0).dot(g1);
133 static_assert(Dune::AlwaysFalse<GEO>::value,
134 "firstDerivative can only be called with wrtCoeff as 0 and 1 indicating first derivative of the "
135 "Green-Lagrange strain w.r.t d and "
136 "alpha, respectively.");
163 template <
int wrtCoeff,
typename GEO,
typename ST,
typename EAST>
164 static auto secondDerivative(
const GEO& geo,
const auto& uFunction,
const auto& localBasis,
const auto& gpIndex,
166 const Eigen::Vector<ST, GEO::mydimension*(GEO::mydimension + 1) / 2>& S,
167 const EAST& easFunction,
const auto& alpha,
const int I = sNaN,
const int J = sNaN) {
168 constexpr int myDim = GEO::mydimension;
169 using namespace Dune::DerivativeDirections;
170 using namespace Dune;
171 static_assert(myDim == 2 or myDim == 3,
172 "An enhancement of displacement gradient can only be performed for the 2D or 3D case.");
173 if constexpr (wrtCoeff == 0) {
174 const auto strainFunction = Dune::greenLagrangeStrains(uFunction);
175 const auto kg = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(I, J)), along(S), on(gridElement));
177 }
else if constexpr (wrtCoeff == 1) {
178 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
179 const typename EAST::HType Harray = easFunction(gpPos);
180 auto kg = Eigen::Matrix<ST, enhancedStrainSize, enhancedStrainSize>::Zero().eval();
182 for (
const auto i : Dune::range(myDim)) {
183 for (
const auto j : Dune::range(myDim)) {
184 const auto sij = 0.5 * S[toVoigt<myDim>(i, j)];
185 for (
const auto p : Dune::range(enhancedStrainSize))
186 for (
const auto q : Dune::range(enhancedStrainSize)) {
187 const auto Hq = Harray[q];
188 kg(p, q) += sij * (Harray[p].col(i).transpose() * Harray[q].col(j) +
189 Harray[p].col(j).transpose() * Harray[q].col(i))
195 }
else if constexpr (wrtCoeff == 2) {
196 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
197 const typename EAST::HType Harray = easFunction(gpPos);
198 auto kg = Eigen::Matrix<ST, enhancedStrainSize, myDim>::Zero().eval();
199 const auto gradUdI = uFunction.evaluateDerivative(gpIndex, wrt(spatialAll, coeff(I)), on(gridElement));
201 std::array<Eigen::Matrix<ST, myDim, myDim>, myDim> dNIdT;
202 std::ranges::fill(dNIdT, Eigen::Matrix<ST, myDim, myDim>::Zero());
203 for (
const auto i : Dune::range(myDim))
204 for (
const auto j : Dune::range(myDim))
205 dNIdT[i](i, j) = getDiagonalEntry(gradUdI[j], 0);
207 for (
const auto i : Dune::range(myDim)) {
208 for (
const auto j : Dune::range(myDim)) {
209 const auto sij = 0.5 * S[toVoigt<myDim>(i, j)];
210 for (
const auto p : Dune::range(enhancedStrainSize))
211 for (
const auto q : Dune::range(myDim)) {
212 kg(p, q) += sij * (Harray[p].col(i).transpose() * dNIdT[q].col(j) +
213 Harray[p].col(j).transpose() * dNIdT[q].col(i))
221 Dune::AlwaysFalse<GEO>::value,
222 "secondDerivative can only be called with wrtCoeff as 0, 1 and 2 indicating second derivative of the "
223 "Green-Lagrange strain w.r.t d and "
224 "alpha and the mixed derivative, respectively.");
228 static constexpr auto name() {
return std::string(
"Displacement Gradient"); }
231 static constexpr int sNaN = std::numeric_limits<int>::signaling_NaN();
233 template <
typename GEO,
typename EAST>
234 static auto computeDisplacementGradient(
const GEO& geo,
const auto& uFunction,
236 const EAST& easFunction,
const auto& alpha) {
237 using ST =
typename std::remove_cvref_t<
decltype(uFunction)>::ctype;
238 constexpr int myDim = GEO::mydimension;
239 using MatrixType = Eigen::Matrix<ST, myDim, myDim>;
241 using namespace Dune::DerivativeDirections;
242 using namespace Dune;
244 const MatrixType Hc = toEigen(uFunction.evaluateDerivative(gpPos, wrt(spatialAll), on(gridElement))).eval();
245 const typename EAST::HType Harray = easFunction(gpPos);
246 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
247 typename EAST::AnsatzType Htilde;
249 for (
const auto p : Dune::range(enhancedStrainSize))
250 Htilde += Harray[p] * alpha[p];
251 const MatrixType H = Hc + Htilde;
255 template <
typename GEO,
typename EAST>
256 static auto computeDeformationGradient(
const GEO& geo,
const auto& uFunction,
258 const EAST& easFunction,
const auto& alpha) {
259 using ST =
typename std::remove_cvref_t<
decltype(uFunction)>::ctype;
260 constexpr int myDim = GEO::mydimension;
261 using MatrixType = Eigen::Matrix<ST, myDim, myDim>;
262 const MatrixType H = computeDisplacementGradient(geo, uFunction, gpPos, easFunction, alpha);
263 const MatrixType F = (MatrixType::Identity() + H).eval();
Helper for the Eigen::Tensor types.
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:182
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: easfunctions/displacementgradient.hh:29
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: easfunctions/displacementgradient.hh:45
static constexpr auto name()
The name of the strain measure enhanced w.r.t EAS method.
Definition: easfunctions/displacementgradient.hh:228
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: easfunctions/displacementgradient.hh:75
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: easfunctions/displacementgradient.hh:164
Definition: utils/dirichletvalues.hh:32