version 0.4
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
29 template <typename Geometry>
30 auto value(const Dune::FieldVector<double, 2>& gpPos, const Geometry& geo, const auto& uFunction) const
31 -> Eigen::Vector3<typename std::remove_cvref_t<decltype(uFunction)>::ctype> {
32 using ScalarType = typename std::remove_cvref_t<decltype(uFunction)>::ctype;
33 Eigen::Vector3<ScalarType> epsV;
34 const auto J = Dune::toEigen(geo.jacobianTransposed(gpPos));
35 using namespace Dune;
36 using namespace Dune::DerivativeDirections;
37 const Eigen::Matrix<ScalarType, 3, 2> gradu = toEigen(
38 uFunction.evaluateDerivative(gpPos, // Here the gpIndex could be passed
39 Dune::wrt(spatialAll), Dune::on(Dune::DerivativeDirections::referenceElement)));
40 const Eigen::Matrix<ScalarType, 2, 3> j = J + gradu.transpose();
41
42 epsV << J.row(0).dot(gradu.col(0)) + 0.5 * gradu.col(0).squaredNorm(),
43 J.row(1).dot(gradu.col(1)) + 0.5 * gradu.col(1).squaredNorm(),
44 j.row(0).dot(j.row(1)) - J.row(0).dot(J.row(1));
45 return epsV;
46 }
47
64 template <typename Geometry, typename ScalarType>
65 auto derivative(const Dune::FieldVector<double, 2>& gpPos, const Eigen::Matrix<ScalarType, 2, 3>& jcur,
66 const auto& dNAtGp, const Geometry& geo, const auto& uFunction, const auto& localBasis,
67 const int node) const {
68 Eigen::Matrix<ScalarType, 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 Geometry, typename ScalarType>
96 auto secondDerivative(const Dune::FieldVector<double, 2>& gpPos, const auto& dNAtGp, const Geometry& geo,
97 const auto& uFunction, const auto& localBasis, const Eigen::Vector3<ScalarType>& S, int I,
98 int J) const {
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 ScalarType NS = dN1i * dN1j * S[0] + dN2i * dN2j * S[1] + (dN1i * dN2j + dN2i * dN1j) * S[2];
104 Eigen::Matrix<ScalarType, 3, 3> kg = Eigen::Matrix<double, 3, 3>::Identity() * NS;
105 return kg;
106 }
107 };
108
109} // namespace Ikarus
Definition: simpleassemblers.hh:21
Definition: resultevaluators.hh:17
Definition: membranestrains.hh:17
auto derivative(const Dune::FieldVector< double, 2 > &gpPos, const Eigen::Matrix< ScalarType, 2, 3 > &jcur, const auto &dNAtGp, const Geometry &geo, const auto &uFunction, const auto &localBasis, const int node) const
Compute the strain-displacement matrix for a given node and integration point.
Definition: membranestrains.hh:65
auto secondDerivative(const Dune::FieldVector< double, 2 > &gpPos, const auto &dNAtGp, const Geometry &geo, const auto &uFunction, const auto &localBasis, const Eigen::Vector3< ScalarType > &S, int I, int J) const
Compute the second derivative of the membrane strain for a given node pair and integration point.
Definition: membranestrains.hh:96
auto value(const Dune::FieldVector< double, 2 > &gpPos, const Geometry &geo, const auto &uFunction) const -> Eigen::Vector3< typename std::remove_cvref_t< decltype(uFunction)>::ctype >
Compute the strain vector at a given integration point.
Definition: membranestrains.hh:30
Definition: resultevaluators.hh:20