11#include <dune/geometry/quadraturerules.hh>
12#include <dune/localfefunctions/cachedlocalBasis/cachedlocalBasis.hh>
13#include <dune/localfefunctions/impl/standardLocalFunction.hh>
25template <
typename PreFE,
typename FE>
26class KirchhoffLoveShell;
36 template <
typename PreFE,
typename FE>
50template <
typename PreFE,
typename FE>
62 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
79 template <
typename ST =
double>
82 Eigen::Matrix<double, 3, 3>
C;
85 Eigen::Matrix<ST, 2, 3>
j;
86 Eigen::Matrix<double, 2, 3>
J;
88 Eigen::Matrix3<double>
H;
89 Eigen::Vector3<ST>
a3N;
90 Eigen::Vector3<ST>
a3;
101 : mat_{pre.material},
102 thickness_{pre.thickness} {}
109 const auto& localView = underlying().localView();
110 assert(localView.bound());
111 const auto& element = localView.element();
112 auto& firstChild = localView.tree().child(0);
114 const auto& fe = firstChild.finiteElement();
116 numberOfNodes_ = fe.size();
117 order_ = fe.localBasis().order();
118 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
119 if (underlying().quadratureRule())
120 localBasis_.bind(underlying().quadratureRule().value(), Dune::bindDerivatives(0, 1, 2));
122 localBasis_.bind(defaultQuadratureRule<myDim>(element, 2 * order_), Dune::bindDerivatives(0, 1, 2));
136 template <
typename ST =
double>
139 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const {
141 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
142 Dune::StandardLocalFunction uFunction(
143 localBasis_, disp, std::make_shared<const Geometry>(underlying().localView().element().
geometry()));
149 [[nodiscard]]
int order()
const {
return order_; }
160 template <
template <
typename,
int,
int>
class RT>
161 requires(supportsResultType<RT>())
165 DUNE_THROW(Dune::NotImplemented,
"No results are implemented");
170 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
171 auto& underlying() {
return static_cast<FE&
>(*this); }
173 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
175 YoungsModulusAndPoissonsRatio mat_;
178 size_t numberOfNodes_{0};
197 const auto& uFunction)
const {
198 using ST =
typename std::remove_cvref_t<
decltype(uFunction)>::ctype;
201 using namespace Dune;
202 using namespace Dune::DerivativeDirections;
203 const auto [X, Jd, Hd] = geo.impl().zeroFirstAndSecondDerivativeOfPosition(gpPos);
206 const Eigen::Matrix<double, 2, 2> A = kin.J * kin.J.transpose();
207 Eigen::Matrix<double, 3, 3> G = Eigen::Matrix<double, 3, 3>::Zero();
209 G.block<2, 2>(0, 0) = A;
211 const Eigen::Matrix<double, 3, 3> GInv = G.inverse();
216 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
218 const auto hessianu = Ndd.transpose().template cast<ST>() * uasMatrix;
219 kin.h = kin.H + hessianu;
220 const Eigen::Matrix<ST, 3, 2> gradu = toEigen(uFunction.evaluateDerivative(
221 gpIndex, Dune::wrt(spatialAll), Dune::on(Dune::DerivativeDirections::referenceElement)));
222 kin.j = kin.J + gradu.transpose();
223 kin.a3N = (kin.j.row(0).cross(kin.j.row(1)));
224 kin.a3 = kin.a3N.normalized();
225 Eigen::Vector<ST, 3> bV = kin.h * kin.a3;
227 const auto BV =
toVoigt(toEigen(geo.impl().secondFundamentalForm(gpPos)));
228 kin.kappaV = BV - bV;
232 template <
typename ST>
235 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const {
237 DUNE_THROW(Dune::NotImplemented,
"MatrixAffordance not implemented: " +
toString(affordance));
238 using namespace Dune::DerivativeDirections;
239 using namespace Dune;
242 const auto geo = underlying().localView().element().geometry();
244 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
245 const auto intElement = geo.integrationElement(gp.position()) * gp.weight();
246 const auto [C, epsV, kappaV, jE, J, h, H, a3N, a3] =
248 const Eigen::Vector<ST, membraneStrainSize> membraneForces = thickness_ * C * epsV;
249 const Eigen::Vector<ST, bendingStrainSize> moments = Dune::power(thickness_, 3) / 12.0 * C * kappaV;
251 const auto& Nd = localBasis_.evaluateJacobian(gpIndex);
252 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
253 for (
size_t i = 0; i < numberOfNodes_; ++i) {
254 Eigen::Matrix<ST, membraneStrainSize, worldDim> bopIMembrane =
257 Eigen::Matrix<ST, bendingStrainSize, worldDim> bopIBending =
bopBending(jE, h, Nd, Ndd, i, a3N, a3);
258 for (
size_t j = i; j < numberOfNodes_; ++j) {
259 auto KBlock = K.template block<worldDim, worldDim>(
worldDim * i,
worldDim * j);
260 Eigen::Matrix<ST, membraneStrainSize, worldDim> bopJMembrane =
262 Eigen::Matrix<ST, bendingStrainSize, worldDim> bopJBending =
bopBending(jE, h, Nd, Ndd, j, a3N, a3);
263 KBlock += thickness_ * bopIMembrane.transpose() * C * bopJMembrane * intElement;
264 KBlock += Dune::power(thickness_, 3) / 12.0 * bopIBending.transpose() * C * bopJBending * intElement;
267 gp.position(), Nd, geo, uFunction, localBasis_, membraneForces, i, j);
268 Eigen::Matrix<ST, worldDim, worldDim> kgBendingIJ =
kgBending(jE, h, Nd, Ndd, a3N, a3, moments, i, j);
269 KBlock += kgMembraneIJ * intElement;
270 KBlock += kgBendingIJ * intElement;
274 K.template triangularView<Eigen::StrictlyLower>() = K.transpose();
277 template <
typename ST>
280 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const {
282 DUNE_THROW(Dune::NotImplemented,
"VectorAffordance not implemented: " +
toString(affordance));
283 using namespace Dune::DerivativeDirections;
284 using namespace Dune;
287 const auto geo = underlying().localView().element().geometry();
290 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
291 const auto [C, epsV, kappaV, jE, J, h, H, a3N, a3] =
293 const Eigen::Vector<ST, 3> membraneForces = thickness_ * C * epsV;
294 const Eigen::Vector<ST, 3> moments = Dune::power(thickness_, 3) / 12.0 * C * kappaV;
296 const auto& Nd = localBasis_.evaluateJacobian(gpIndex);
297 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
298 for (
size_t i = 0; i < numberOfNodes_; ++i) {
299 Eigen::Matrix<ST, 3, 3> bopIMembrane =
301 Eigen::Matrix<ST, 3, 3> bopIBending =
bopBending(jE, h, Nd, Ndd, i, a3N, a3);
302 force.template segment<3>(3 * i) +=
303 bopIMembrane.transpose() * membraneForces * geo.integrationElement(gp.position()) * gp.weight();
304 force.template segment<3>(3 * i) +=
305 bopIBending.transpose() * moments * geo.integrationElement(gp.position()) * gp.weight();
310 template <
typename ST>
313 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const -> ST {
315 DUNE_THROW(Dune::NotImplemented,
"ScalarAffordance not implemented: " +
toString(affordance));
316 using namespace Dune::DerivativeDirections;
317 using namespace Dune;
319 const auto& lambda = par.parameter();
322 const auto geo = underlying().localView().element().geometry();
324 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
325 const auto [C, epsV, kappaV, j, J, h, H, a3N, a3] =
328 const ST membraneEnergy = 0.5 * thickness_ * epsV.dot(C * epsV);
329 const ST bendingEnergy = 0.5 * Dune::power(thickness_, 3) / 12.0 * kappaV.dot(C * kappaV);
330 energy += (membraneEnergy + bendingEnergy) * geo.integrationElement(gp.position()) * gp.weight();
336 template <
typename ST>
337 Eigen::Matrix<ST, 3, 3>
kgBending(
const Eigen::Matrix<ST, 2, 3>& jcur,
const Eigen::Matrix3<ST>& h,
const auto& dN,
338 const auto& ddN,
const Eigen::Vector3<ST>& a3N,
const Eigen::Vector3<ST>& a3,
339 const Eigen::Vector3<ST>& S,
int I,
int J)
const {
340 Eigen::Matrix<ST, 3, 3> kg;
343 const auto& dN1i = dN(I, 0);
344 const auto& dN1j = dN(J, 0);
345 const auto& dN2i = dN(I, 1);
346 const auto& dN2j = dN(J, 1);
348 const Eigen::Matrix<ST, 3, 3> P =
349 1.0 / a3N.norm() * (Eigen::Matrix<double, 3, 3>::Identity() - a3 * a3.transpose());
351 const auto a1dxI = Eigen::Matrix<double, 3, 3>::Identity() *
354 const auto a2dxI = Eigen::Matrix<double, 3, 3>::Identity() * dN2i;
355 const auto a1dxJ = Eigen::Matrix<double, 3, 3>::Identity() * dN1j;
356 const auto a2dxJ = Eigen::Matrix<double, 3, 3>::Identity() * dN2j;
357 const auto a1 = jcur.row(0);
358 const auto a2 = jcur.row(1);
359 const Eigen::Matrix<ST, 3, 3> a3NdI = a1dxI.colwise().cross(a2) - a2dxI.colwise().cross(a1);
360 const Eigen::Matrix<ST, 3, 3> a3NdJ = a1dxJ.colwise().cross(a2) - a2dxJ.colwise().cross(a1);
361 Eigen::Matrix<ST, 3, 3> a3dI = P * a3NdI;
362 Eigen::Matrix<ST, 3, 3> a3dJ = P * a3NdJ;
363 for (
int i = 0; i < 3; ++i) {
364 const auto a_albe = h.row(i).transpose();
365 const auto& ddNI = ddN(I, i);
366 const auto& ddNJ = ddN(J, i);
367 Eigen::Vector3<ST> vecd = P * a_albe;
369 Eigen::Matrix<ST, 3, 3> a3Ndd =
370 1.0 / a3N.squaredNorm() *
371 ((3 * a3 * a3.transpose() - Eigen::Matrix<double, 3, 3>::Identity()) * (a3.dot(a_albe)) -
372 a_albe * a3.transpose() - a3 * a_albe.transpose());
374 Eigen::Matrix<ST, 3, 3> secondDerivativeDirectorIJ =
skew(((dN2i * dN1j - dN1i * dN2j) * vecd).eval());
375 kg -= (a3NdI.transpose() * a3Ndd * a3NdJ + secondDerivativeDirectorIJ + (ddNI * a3dJ + ddNJ * a3dI.transpose())) *
376 S[i] * (i == 2 ? 2 : 1);
382 template <
typename ST>
383 Eigen::Matrix<ST, 3, 3>
bopBending(
const Eigen::Matrix<ST, 2, 3>& jcur,
const Eigen::Matrix3<ST>& h,
const auto& dN,
384 const auto& ddN,
const int node,
const Eigen::Vector3<ST>& a3N,
385 const Eigen::Vector3<ST>& a3)
const {
386 const Eigen::Matrix<ST, 3, 3> a1dxI =
387 Eigen::Matrix<double, 3, 3>::Identity() * dN(node, 0);
389 const Eigen::Matrix<ST, 3, 3> a2dxI = Eigen::Matrix<double, 3, 3>::Identity() * dN(node, 1);
390 const auto a1 = jcur.row(0);
391 const auto a2 = jcur.row(1);
392 const Eigen::Matrix<ST, 3, 3> a3NdI =
393 a1dxI.colwise().cross(a2) - a2dxI.colwise().cross(a1);
395 const Eigen::Matrix<ST, 3, 3> a3d1 =
396 1.0 / a3N.norm() * (Eigen::Matrix<double, 3, 3>::Identity() - a3 * a3.transpose()) * a3NdI;
398 Eigen::Matrix<ST, 3, 3> bop = -(h * a3d1 + (a3 * ddN.row(node)).
transpose());
410 Eigen::Matrix<double, 3, 3>
materialTangent(
const Eigen::Matrix<double, 3, 3>& Aconv)
const {
411 const double lambda = mat_.emodul * mat_.nu / ((1.0 + mat_.nu) * (1.0 - 2.0 * mat_.nu));
412 const double mu = mat_.emodul / (2.0 * (1.0 + mat_.nu));
413 const double lambdbar = 2.0 * lambda * mu / (lambda + 2.0 * mu);
414 Eigen::TensorFixedSize<double, Eigen::Sizes<3, 3, 3, 3>> moduli;
415 const auto AconvT =
tensorView(Aconv, std::array<Eigen::Index, 2>({3, 3}));
416 moduli = lambdbar *
dyadic(AconvT, AconvT).eval() + 2.0 * mu * symmetricFourthOrder<double>(Aconv, Aconv);
419 Eigen::Matrix<double, 3, 3> C33 = C({0, 1, 5}, {0, 1, 5});
Free function for creating a tensor product quadrature rule and other Dune::Quadrature rule utils.
Helper for the autodiff library.
Implementation of membrane strain for shells.
Header file for types of loads in Ikarus finite element mechanics.
Definitions of ResultTypes used for finite element results.
Material property functions and conversion utilities.
Definition of the LinearElastic class for finite element mechanics computations.
Derived skew(const Eigen::MatrixBase< Derived > &A)
Returns the skew part of a matrix.
Definition: linearalgebrahelper.hh:434
auto viewAsEigenMatrixAsDynFixed(Dune::BlockVector< ValueType > &blockedVector)
View Dune::BlockVector as an Eigen::Matrix with dynamic rows and fixed columns depending on the size ...
Definition: linearalgebrahelper.hh:87
Eigen::Tensor< typename Derived::Scalar, rank > tensorView(const Eigen::EigenBase< Derived > &matrix, const std::array< T, rank > &dims)
View an Eigen matrix as an Eigen Tensor with specified dimensions.
Definition: tensorutils.hh:32
auto dyadic(const auto &A_ij, const auto &B_kl)
Computes the dyadic product of two Eigen tensors.
Definition: tensorutils.hh:47
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:182
Definition: assemblermanipulatorbuildingblocks.hh:22
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:65
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:50
auto transpose(const Eigen::EigenBase< Derived > &A)
auto kirchhoffLoveShell(const KlArgs &args)
A helper function to create a Kirchhoff-Love shell pre finite element.
Definition: kirchhoffloveshell.hh:440
constexpr std::string toString(DBCOption _e)
Definition: dirichletbcenforcement.hh:8
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:39
@ mechanicalPotentialEnergy
Definition: utils/dirichletvalues.hh:35
FE class is a base class for all finite elements.
Definition: febase.hh:79
FETraits< BH, useEigenRef, useFlat > Traits
Definition: febase.hh:38
Class representing the requirements for finite element calculations.
Definition: ferequirements.hh:224
const SolutionVectorType & globalSolution() const
Get the global solution vector.
Definition: ferequirements.hh:314
const PM & parameter() const
Get the parameter value.
Definition: ferequirements.hh:349
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:164
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:277
Traits for handling finite elements.
Definition: fetraits.hh:25
typename Basis::LocalView LocalView
Type of the local view.
Definition: fetraits.hh:42
typename Element::Geometry Geometry
Type of the element geometry.
Definition: fetraits.hh:51
BH BasisHandler
Type of the basis of the finite element.
Definition: fetraits.hh:27
typename Basis::GridView GridView
Type of the grid view.
Definition: fetraits.hh:45
static constexpr int worlddim
Dimension of the world space.
Definition: fetraits.hh:60
typename BasisHandler::FlatBasis FlatBasis
Type of the flat basis.
Definition: fetraits.hh:33
typename LocalView::Element Element
Type of the grid element.
Definition: fetraits.hh:48
static constexpr int mydim
Dimension of the geometry.
Definition: fetraits.hh:63
Kirchhoff-Love shell finite element class.
Definition: kirchhoffloveshell.hh:52
Eigen::Matrix< double, 3, 3 > materialTangent(const Eigen::Matrix< double, 3, 3 > &Aconv) const
Gets the material tangent matrix for the linear elastic material.
Definition: kirchhoffloveshell.hh:410
typename Traits::Element Element
Definition: kirchhoffloveshell.hh:61
typename Traits::Geometry Geometry
Definition: kirchhoffloveshell.hh:59
void calculateVectorImpl(const Requirement &par, const VectorAffordance &affordance, typename Traits::template VectorType< ST > force, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const
Definition: kirchhoffloveshell.hh:278
auto displacementFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const
Get the displacement function and nodal displacements.
Definition: kirchhoffloveshell.hh:137
auto computeMaterialAndStrains(const Dune::FieldVector< double, 2 > &gpPos, int gpIndex, const Geometry &geo, const auto &uFunction) const
Compute material properties and strains at a given integration point.
Definition: kirchhoffloveshell.hh:196
size_t numberOfNodes() const
Definition: kirchhoffloveshell.hh:148
static constexpr int myDim
Definition: kirchhoffloveshell.hh:65
static constexpr int bendingStrainSize
Definition: kirchhoffloveshell.hh:68
auto calculateAtImpl(const Requirement &req, const Dune::FieldVector< double, Traits::mydim > &local) -> ResultWrapper< RT< double, myDim, worldDim >, ResultShape::Vector >
Calculates a requested result at a specific local position.
Definition: kirchhoffloveshell.hh:162
Geometry geometry() const
Definition: kirchhoffloveshell.hh:147
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: kirchhoffloveshell.hh:62
typename Traits::GridView GridView
Definition: kirchhoffloveshell.hh:60
static constexpr int worldDim
Definition: kirchhoffloveshell.hh:66
void bindImpl()
A helper function to bind the local view to the element.
Definition: kirchhoffloveshell.hh:108
typename Traits::FlatBasis FlatBasis
Definition: kirchhoffloveshell.hh:56
KirchhoffLoveShell(const Pre &pre)
Constructor for the KirchhoffLoveShell class.
Definition: kirchhoffloveshell.hh:100
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType< ST > K, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const
Definition: kirchhoffloveshell.hh:233
typename Traits::LocalView LocalView
Definition: kirchhoffloveshell.hh:58
typename Traits::BasisHandler BasisHandler
Definition: kirchhoffloveshell.hh:55
Eigen::Matrix< ST, 3, 3 > kgBending(const Eigen::Matrix< ST, 2, 3 > &jcur, const Eigen::Matrix3< ST > &h, const auto &dN, const auto &ddN, const Eigen::Vector3< ST > &a3N, const Eigen::Vector3< ST > &a3, const Eigen::Vector3< ST > &S, int I, int J) const
Definition: kirchhoffloveshell.hh:337
Eigen::Matrix< ST, 3, 3 > bopBending(const Eigen::Matrix< ST, 2, 3 > &jcur, const Eigen::Matrix3< ST > &h, const auto &dN, const auto &ddN, const int node, const Eigen::Vector3< ST > &a3N, const Eigen::Vector3< ST > &a3) const
Definition: kirchhoffloveshell.hh:383
static constexpr int membraneStrainSize
Definition: kirchhoffloveshell.hh:67
auto calculateScalarImpl(const Requirement &par, const ScalarAffordance &affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const -> ST
Definition: kirchhoffloveshell.hh:311
int order() const
Definition: kirchhoffloveshell.hh:149
A PreFE struct for Kirchhoff-Love shell elements.
Definition: kirchhoffloveshell.hh:32
double thickness
Definition: kirchhoffloveshell.hh:34
YoungsModulusAndPoissonsRatio material
Definition: kirchhoffloveshell.hh:33
A structure representing kinematic variables.
Definition: kirchhoffloveshell.hh:81
Eigen::Matrix3< double > H
Hessian of the reference geometry.
Definition: kirchhoffloveshell.hh:88
Eigen::Vector3< ST > kappaV
bending strain in Voigt notation
Definition: kirchhoffloveshell.hh:84
Eigen::Matrix< double, 3, 3 > C
material tangent
Definition: kirchhoffloveshell.hh:82
Eigen::Matrix< ST, 2, 3 > j
Jacobian of the deformed geometry.
Definition: kirchhoffloveshell.hh:85
Eigen::Vector3< ST > epsV
membrane strain in Voigt notation
Definition: kirchhoffloveshell.hh:83
Eigen::Matrix3< ST > h
Hessian of the deformed geometry.
Definition: kirchhoffloveshell.hh:87
Eigen::Vector3< ST > a3N
Normal vector of the deformed geometry.
Definition: kirchhoffloveshell.hh:89
Eigen::Vector3< ST > a3
normalized normal vector of the deformed geometry
Definition: kirchhoffloveshell.hh:90
Eigen::Matrix< double, 2, 3 > J
Jacobian of the reference geometry.
Definition: kirchhoffloveshell.hh:86
A struct containing information about the Youngs Modulus, Poisson's ratio and the thickness for the K...
Definition: kirchhoffloveshell.hh:429
double youngs_modulus
Definition: kirchhoffloveshell.hh:430
double nu
Definition: kirchhoffloveshell.hh:431
double thickness
Definition: kirchhoffloveshell.hh:432
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:64
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:95
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:30
see https://en.wikipedia.org/wiki/Lam%C3%A9_parameters Structure representing Young's modulus and she...
Definition: physicshelper.hh:19
Definition: utils/dirichletvalues.hh:37