version 0.4.1
membranestrains.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10#include <ranges>
11
12#include <dune/common/fvector.hh>
13
14#include <Eigen/Core>
15namespace Ikarus {
16
18{
30 template <typename GEO>
31 static auto value(const Dune::FieldVector<double, 2>& gpPos, const GEO& geo,
32 const auto& uFunction) -> Eigen::Vector3<typename std::remove_cvref_t<decltype(uFunction)>::ctype> {
33 using ScalarType = typename std::remove_cvref_t<decltype(uFunction)>::ctype;
34 Eigen::Vector3<ScalarType> epsV;
35 const auto J = Dune::toEigen(geo.jacobianTransposed(gpPos));
36 using namespace Dune;
37 using namespace Dune::DerivativeDirections;
38 const Eigen::Matrix<ScalarType, 3, 2> gradu = toEigen(
39 uFunction.evaluateDerivative(gpPos, // Here the gpIndex could be passed
40 Dune::wrt(spatialAll), Dune::on(Dune::DerivativeDirections::referenceElement)));
41 const Eigen::Matrix<ScalarType, 2, 3> j = J + gradu.transpose();
42
43 epsV << J.row(0).dot(gradu.col(0)) + 0.5 * gradu.col(0).squaredNorm(),
44 J.row(1).dot(gradu.col(1)) + 0.5 * gradu.col(1).squaredNorm(), j.row(0).dot(j.row(1)) - J.row(0).dot(J.row(1));
45 return epsV;
46 }
47
64 template <typename GEO, typename ST>
65 static auto derivative(const Dune::FieldVector<double, 2>& gpPos, const Eigen::Matrix<ST, 2, 3>& jcur,
66 const auto& dNAtGp, const GEO& geo, const auto& uFunction, const auto& localBasis,
67 const int node) {
68 Eigen::Matrix<ST, 3, 3> bop;
69 bop.row(0) = jcur.row(0) * dNAtGp(node, 0);
70 bop.row(1) = jcur.row(1) * dNAtGp(node, 1);
71 bop.row(2) = jcur.row(0) * dNAtGp(node, 1) + jcur.row(1) * dNAtGp(node, 0);
72
73 return bop;
74 }
75
95 template <typename GEO, typename ST>
96 static auto secondDerivative(const Dune::FieldVector<double, 2>& gpPos, const auto& dNAtGp, const GEO& geo,
97 const auto& uFunction, const auto& localBasis, const Eigen::Vector3<ST>& S, int I,
98 int J) {
99 const auto& dN1i = dNAtGp(I, 0);
100 const auto& dN1j = dNAtGp(J, 0);
101 const auto& dN2i = dNAtGp(I, 1);
102 const auto& dN2j = dNAtGp(J, 1);
103 const ST NS = dN1i * dN1j * S[0] + dN2i * dN2j * S[1] + (dN1i * dN2j + dN2i * dN1j) * S[2];
104 Eigen::Matrix<ST, 3, 3> kg = Eigen::Matrix<double, 3, 3>::Identity() * NS;
105 return kg;
106 }
107};
108
109} // namespace Ikarus
Definition: simpleassemblers.hh:22
Definition: utils/dirichletvalues.hh:28
Definition: membranestrains.hh:18
static auto derivative(const Dune::FieldVector< double, 2 > &gpPos, const Eigen::Matrix< ST, 2, 3 > &jcur, const auto &dNAtGp, const GEO &geo, const auto &uFunction, const auto &localBasis, const int node)
Compute the strain-displacement matrix for a given node and integration point.
Definition: membranestrains.hh:65
static auto secondDerivative(const Dune::FieldVector< double, 2 > &gpPos, const auto &dNAtGp, const GEO &geo, const auto &uFunction, const auto &localBasis, const Eigen::Vector3< ST > &S, int I, int J)
Compute the second derivative of the membrane strain for a given node pair and integration point.
Definition: membranestrains.hh:96
static auto value(const Dune::FieldVector< double, 2 > &gpPos, const GEO &geo, const auto &uFunction) -> Eigen::Vector3< typename std::remove_cvref_t< decltype(uFunction)>::ctype >
Compute the strain vector at a given integration point.
Definition: membranestrains.hh:31
Definition: utils/dirichletvalues.hh:30