version 0.4.1
kirchhoffloveshell.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
11#include <dune/geometry/quadraturerules.hh>
12#include <dune/localfefunctions/cachedlocalBasis/cachedlocalBasis.hh>
13#include <dune/localfefunctions/impl/standardLocalFunction.hh>
14
22
23namespace Ikarus {
24template <typename PreFE, typename FE>
25class KirchhoffLoveShell;
26
31{
33 double thickness;
34
35 template <typename PreFE, typename FE>
37};
38
49template <typename PreFE, typename FE>
51{
52public:
55 using FlatBasis = typename Traits::FlatBasis;
58 using LocalView = typename Traits::LocalView;
59 using Geometry = typename Traits::Geometry;
60 using GridView = typename Traits::GridView;
61 using Element = typename Traits::Element;
62 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
64
65 static constexpr int myDim = Traits::mydim;
66 static constexpr int worldDim = Traits::worlddim;
67 static constexpr int membraneStrainSize = 3;
68 static constexpr int bendingStrainSize = 3;
69
79 template <typename ST = double>
81 {
82 Eigen::Matrix<double, 3, 3> C;
83 Eigen::Vector3<ST> epsV;
84 Eigen::Vector3<ST> kappaV;
85 Eigen::Matrix<ST, 2, 3> j;
86 Eigen::Matrix<double, 2, 3> J;
87 Eigen::Matrix3<ST> h;
88 Eigen::Matrix3<double> H;
89 Eigen::Vector3<ST> a3N;
90 Eigen::Vector3<ST> a3;
91 };
92
101 : mat_{pre.material},
102 thickness_{pre.thickness} {}
103
104protected:
108 void bindImpl() {
109 const auto& localView = underlying().localView();
110 assert(localView.bound());
111 const auto& element = localView.element();
112 auto& firstChild = localView.tree().child(0);
113
114 const auto& fe = firstChild.finiteElement();
115 // geo_.emplace(element.geometry());
116 numberOfNodes_ = fe.size();
117 order_ = 2 * (fe.localBasis().order());
118 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
119 if constexpr (requires { element.impl().getQuadratureRule(order_); })
120 if (element.impl().isTrimmed())
121 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1, 2));
122 else
123 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
124 Dune::bindDerivatives(0, 1, 2));
125 else
126 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
127 Dune::bindDerivatives(0, 1, 2));
128 }
129
130public:
141 template <typename ST = double>
143 const Requirement& par,
144 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const {
145 const auto& d = par.globalSolution();
146 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
147 Dune::StandardLocalFunction uFunction(
148 localBasis_, disp, std::make_shared<const Geometry>(underlying().localView().element().geometry()));
149 return uFunction;
150 }
151
152 Geometry geometry() const { return underlying().localView().element().geometry(); }
153 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
154 [[nodiscard]] int order() const { return order_; }
155
165 template <template <typename, int, int> class RT>
166 requires(supportsResultType<RT>())
167 auto calculateAtImpl([[maybe_unused]] const Requirement& req,
168 [[maybe_unused]] const Dune::FieldVector<double, Traits::mydim>& local)
170 DUNE_THROW(Dune::NotImplemented, "No results are implemented");
171 }
172
173private:
174 //> CRTP
175 const auto& underlying() const { return static_cast<const FE&>(*this); }
176 auto& underlying() { return static_cast<FE&>(*this); }
177 // std::optional<Geometry> geo_;
178 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
179 // DefaultMembraneStrain membraneStrain_;
180 YoungsModulusAndPoissonsRatio mat_;
181 double thickness_;
182
183 size_t numberOfNodes_{0};
184 int order_{};
185
186protected:
201 auto computeMaterialAndStrains(const Dune::FieldVector<double, 2>& gpPos, int gpIndex, const Geometry& geo,
202 const auto& uFunction) const {
203 using ST = typename std::remove_cvref_t<decltype(uFunction)>::ctype;
204
206 using namespace Dune;
207 using namespace Dune::DerivativeDirections;
208 const auto [X, Jd, Hd] = geo.impl().zeroFirstAndSecondDerivativeOfPosition(gpPos);
209 kin.J = toEigen(Jd);
210 kin.H = toEigen(Hd);
211 const Eigen::Matrix<double, 2, 2> A = kin.J * kin.J.transpose();
212 Eigen::Matrix<double, 3, 3> G = Eigen::Matrix<double, 3, 3>::Zero();
213
214 G.block<2, 2>(0, 0) = A;
215 G(2, 2) = 1;
216 const Eigen::Matrix<double, 3, 3> GInv = G.inverse();
217 kin.C = materialTangent(GInv);
218
219 kin.epsV = DefaultMembraneStrain::value(gpPos, geo, uFunction);
220
221 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
222 const auto uasMatrix = Dune::viewAsEigenMatrixAsDynFixed(uFunction.coefficientsRef());
223 const auto hessianu = Ndd.transpose().template cast<ST>() * uasMatrix;
224 kin.h = kin.H + hessianu;
225 const Eigen::Matrix<ST, 3, 2> gradu = toEigen(uFunction.evaluateDerivative(
226 gpIndex, Dune::wrt(spatialAll), Dune::on(Dune::DerivativeDirections::referenceElement)));
227 kin.j = kin.J + gradu.transpose();
228 kin.a3N = (kin.j.row(0).cross(kin.j.row(1)));
229 kin.a3 = kin.a3N.normalized();
230 Eigen::Vector<ST, 3> bV = kin.h * kin.a3;
231 bV(2) *= 2; // Voigt notation requires the two here
232 const auto BV = toVoigt(toEigen(geo.impl().secondFundamentalForm(gpPos)));
233 kin.kappaV = BV - bV;
234 return kin;
235 }
236
237 template <typename ST>
239 const Requirement& par, const MatrixAffordance& affordance, typename Traits::template MatrixType<ST> K,
240 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const {
241 if (affordance != MatrixAffordance::stiffness)
242 DUNE_THROW(Dune::NotImplemented, "MatrixAffordance not implemented: " + toString(affordance));
243 using namespace Dune::DerivativeDirections;
244 using namespace Dune;
245 const auto uFunction = displacementFunction(par, dx);
246 const auto& lambda = par.parameter();
247 const auto geo = underlying().localView().element().geometry();
248
249 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
250 const auto intElement = geo.integrationElement(gp.position()) * gp.weight();
251 const auto [C, epsV, kappaV, jE, J, h, H, a3N, a3] =
252 computeMaterialAndStrains(gp.position(), gpIndex, geo, uFunction);
253 const Eigen::Vector<ST, membraneStrainSize> membraneForces = thickness_ * C * epsV;
254 const Eigen::Vector<ST, bendingStrainSize> moments = Dune::power(thickness_, 3) / 12.0 * C * kappaV;
255
256 const auto& Nd = localBasis_.evaluateJacobian(gpIndex);
257 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
258 for (size_t i = 0; i < numberOfNodes_; ++i) {
259 Eigen::Matrix<ST, membraneStrainSize, worldDim> bopIMembrane =
260 DefaultMembraneStrain::derivative(gp.position(), jE, Nd, geo, uFunction, localBasis_, i);
261
262 Eigen::Matrix<ST, bendingStrainSize, worldDim> bopIBending = bopBending(jE, h, Nd, Ndd, i, a3N, a3);
263 for (size_t j = i; j < numberOfNodes_; ++j) {
264 auto KBlock = K.template block<worldDim, worldDim>(worldDim * i, worldDim * j);
265 Eigen::Matrix<ST, membraneStrainSize, worldDim> bopJMembrane =
266 DefaultMembraneStrain::derivative(gp.position(), jE, Nd, geo, uFunction, localBasis_, j);
267 Eigen::Matrix<ST, bendingStrainSize, worldDim> bopJBending = bopBending(jE, h, Nd, Ndd, j, a3N, a3);
268 KBlock += thickness_ * bopIMembrane.transpose() * C * bopJMembrane * intElement;
269 KBlock += Dune::power(thickness_, 3) / 12.0 * bopIBending.transpose() * C * bopJBending * intElement;
270
271 Eigen::Matrix<ST, worldDim, worldDim> kgMembraneIJ = DefaultMembraneStrain::secondDerivative(
272 gp.position(), Nd, geo, uFunction, localBasis_, membraneForces, i, j);
273 Eigen::Matrix<ST, worldDim, worldDim> kgBendingIJ = kgBending(jE, h, Nd, Ndd, a3N, a3, moments, i, j);
274 KBlock += kgMembraneIJ * intElement;
275 KBlock += kgBendingIJ * intElement;
276 }
277 }
278 }
279 K.template triangularView<Eigen::StrictlyLower>() = K.transpose();
280 }
281
282 template <typename ST>
284 const Requirement& par, const VectorAffordance& affordance, typename Traits::template VectorType<ST> force,
285 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const {
286 if (affordance != VectorAffordance::forces)
287 DUNE_THROW(Dune::NotImplemented, "VectorAffordance not implemented: " + toString(affordance));
288 using namespace Dune::DerivativeDirections;
289 using namespace Dune;
290 const auto uFunction = displacementFunction(par, dx);
291 const auto& lambda = par.parameter();
292 const auto geo = underlying().localView().element().geometry();
293
294 // Internal forces
295 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
296 const auto [C, epsV, kappaV, jE, J, h, H, a3N, a3] =
297 computeMaterialAndStrains(gp.position(), gpIndex, geo, uFunction);
298 const Eigen::Vector<ST, 3> membraneForces = thickness_ * C * epsV;
299 const Eigen::Vector<ST, 3> moments = Dune::power(thickness_, 3) / 12.0 * C * kappaV;
300
301 const auto& Nd = localBasis_.evaluateJacobian(gpIndex);
302 const auto& Ndd = localBasis_.evaluateSecondDerivatives(gpIndex);
303 for (size_t i = 0; i < numberOfNodes_; ++i) {
304 Eigen::Matrix<ST, 3, 3> bopIMembrane =
305 DefaultMembraneStrain::derivative(gp.position(), jE, Nd, geo, uFunction, localBasis_, i);
306 Eigen::Matrix<ST, 3, 3> bopIBending = bopBending(jE, h, Nd, Ndd, i, a3N, a3);
307 force.template segment<3>(3 * i) +=
308 bopIMembrane.transpose() * membraneForces * geo.integrationElement(gp.position()) * gp.weight();
309 force.template segment<3>(3 * i) +=
310 bopIBending.transpose() * moments * geo.integrationElement(gp.position()) * gp.weight();
311 }
312 }
313 }
314
315 template <typename ST>
317 const Requirement& par, const ScalarAffordance& affordance,
318 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const -> ST {
320 DUNE_THROW(Dune::NotImplemented, "ScalarAffordance not implemented: " + toString(affordance));
321 using namespace Dune::DerivativeDirections;
322 using namespace Dune;
323 const auto uFunction = displacementFunction(par, dx);
324 const auto& lambda = par.parameter();
325 ST energy = 0.0;
326
327 const auto geo = underlying().localView().element().geometry();
328
329 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
330 const auto [C, epsV, kappaV, j, J, h, H, a3N, a3] =
331 computeMaterialAndStrains(gp.position(), gpIndex, geo, uFunction);
332
333 const ST membraneEnergy = 0.5 * thickness_ * epsV.dot(C * epsV);
334 const ST bendingEnergy = 0.5 * Dune::power(thickness_, 3) / 12.0 * kappaV.dot(C * kappaV);
335 energy += (membraneEnergy + bendingEnergy) * geo.integrationElement(gp.position()) * gp.weight();
336 }
337
338 return energy;
339 }
340
341 template <typename ST>
342 Eigen::Matrix<ST, 3, 3> kgBending(const Eigen::Matrix<ST, 2, 3>& jcur, const Eigen::Matrix3<ST>& h, const auto& dN,
343 const auto& ddN, const Eigen::Vector3<ST>& a3N, const Eigen::Vector3<ST>& a3,
344 const Eigen::Vector3<ST>& S, int I, int J) const {
345 Eigen::Matrix<ST, 3, 3> kg;
346 kg.setZero();
347
348 const auto& dN1i = dN(I, 0);
349 const auto& dN1j = dN(J, 0);
350 const auto& dN2i = dN(I, 1);
351 const auto& dN2j = dN(J, 1);
352
353 const Eigen::Matrix<ST, 3, 3> P =
354 1.0 / a3N.norm() * (Eigen::Matrix<double, 3, 3>::Identity() - a3 * a3.transpose());
355
356 const auto a1dxI = Eigen::Matrix<double, 3, 3>::Identity() *
357 dN1i; // the auto here allows the exploitation of the identity matrices,
358 // due to Eigen's expression templates
359 const auto a2dxI = Eigen::Matrix<double, 3, 3>::Identity() * dN2i;
360 const auto a1dxJ = Eigen::Matrix<double, 3, 3>::Identity() * dN1j;
361 const auto a2dxJ = Eigen::Matrix<double, 3, 3>::Identity() * dN2j;
362 const auto a1 = jcur.row(0);
363 const auto a2 = jcur.row(1);
364 const Eigen::Matrix<ST, 3, 3> a3NdI = a1dxI.colwise().cross(a2) - a2dxI.colwise().cross(a1);
365 const Eigen::Matrix<ST, 3, 3> a3NdJ = a1dxJ.colwise().cross(a2) - a2dxJ.colwise().cross(a1);
366 Eigen::Matrix<ST, 3, 3> a3dI = P * a3NdI;
367 Eigen::Matrix<ST, 3, 3> a3dJ = P * a3NdJ;
368 for (int i = 0; i < 3; ++i) {
369 const auto a_albe = h.row(i).transpose();
370 const auto& ddNI = ddN(I, i);
371 const auto& ddNJ = ddN(J, i);
372 Eigen::Vector3<ST> vecd = P * a_albe;
373
374 Eigen::Matrix<ST, 3, 3> a3Ndd =
375 1.0 / a3N.squaredNorm() *
376 ((3 * a3 * a3.transpose() - Eigen::Matrix<double, 3, 3>::Identity()) * (a3.dot(a_albe)) -
377 a_albe * a3.transpose() - a3 * a_albe.transpose());
378
379 Eigen::Matrix<ST, 3, 3> secondDerivativeDirectorIJ = skew(((dN2i * dN1j - dN1i * dN2j) * vecd).eval());
380 kg -= (a3NdI.transpose() * a3Ndd * a3NdJ + secondDerivativeDirectorIJ + (ddNI * a3dJ + ddNJ * a3dI.transpose())) *
381 S[i] * (i == 2 ? 2 : 1);
382 }
383
384 return kg;
385 }
386
387 template <typename ST>
388 Eigen::Matrix<ST, 3, 3> bopBending(const Eigen::Matrix<ST, 2, 3>& jcur, const Eigen::Matrix3<ST>& h, const auto& dN,
389 const auto& ddN, const int node, const Eigen::Vector3<ST>& a3N,
390 const Eigen::Vector3<ST>& a3) const {
391 const Eigen::Matrix<ST, 3, 3> a1dxI =
392 Eigen::Matrix<double, 3, 3>::Identity() * dN(node, 0); // this should be double
393 // but the cross-product below complains otherwise
394 const Eigen::Matrix<ST, 3, 3> a2dxI = Eigen::Matrix<double, 3, 3>::Identity() * dN(node, 1);
395 const auto a1 = jcur.row(0);
396 const auto a2 = jcur.row(1);
397 const Eigen::Matrix<ST, 3, 3> a3NdI =
398 a1dxI.colwise().cross(a2) - a2dxI.colwise().cross(a1); // minus needed since order has
399 // to be swapped to get column-wise cross product working
400 const Eigen::Matrix<ST, 3, 3> a3d1 =
401 1.0 / a3N.norm() * (Eigen::Matrix<double, 3, 3>::Identity() - a3 * a3.transpose()) * a3NdI;
402
403 Eigen::Matrix<ST, 3, 3> bop = -(h * a3d1 + (a3 * ddN.row(node)).transpose());
404 bop.row(2) *= 2;
405
406 return bop;
407 }
408
415 Eigen::Matrix<double, 3, 3> materialTangent(const Eigen::Matrix<double, 3, 3>& Aconv) const {
416 const double lambda = mat_.emodul * mat_.nu / ((1.0 + mat_.nu) * (1.0 - 2.0 * mat_.nu));
417 const double mu = mat_.emodul / (2.0 * (1.0 + mat_.nu));
418 const double lambdbar = 2.0 * lambda * mu / (lambda + 2.0 * mu);
419 Eigen::TensorFixedSize<double, Eigen::Sizes<3, 3, 3, 3>> moduli;
420 const auto AconvT = tensorView(Aconv, std::array<Eigen::Index, 2>({3, 3}));
421 moduli = lambdbar * dyadic(AconvT, AconvT).eval() + 2.0 * mu * symmetricFourthOrder<double>(Aconv, Aconv);
422
423 auto C = toVoigt(moduli);
424 Eigen::Matrix<double, 3, 3> C33 = C({0, 1, 5}, {0, 1, 5});
425 return C33;
426 }
427};
428
433struct KlArgs
434{
436 double nu;
437 double thickness;
438};
439
445auto kirchhoffLoveShell(const KlArgs& args) {
446 KirchhoffLoveShellPre pre({.emodul = args.youngs_modulus, .nu = args.nu}, args.thickness);
447
448 return pre;
449}
450} // namespace Ikarus
Helper for the autodiff library.
Definitions of ResultTypes used for finite element results.
Material property functions and conversion utilities.
Definition of the LinearElastic class for finite element mechanics computations.
Implementation of membrane strain for shells.
Header file for types of loads in Ikarus finite element mechanics.
Derived skew(const Eigen::MatrixBase< Derived > &A)
Returns the skew part of a matrix.
Definition: linearalgebrahelper.hh:410
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:88
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:166
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
Definition: assemblermanipulatorbuildingblocks.hh:22
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
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:445
constexpr std::string toString(DBCOption _e)
Definition: dirichletbcenforcement.hh:7
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
Definition: utils/dirichletvalues.hh:30
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:252
SolutionVectorReturnType globalSolution()
Get the global solution vector.
Definition: ferequirements.hh:308
PMHelper::ConstReturnType parameter() const
Get the parameter value.
Definition: ferequirements.hh:325
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:159
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:272
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:51
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:415
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:283
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:142
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:201
size_t numberOfNodes() const
Definition: kirchhoffloveshell.hh:153
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:167
Geometry geometry() const
Definition: kirchhoffloveshell.hh:152
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:55
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:238
typename Traits::LocalView LocalView
Definition: kirchhoffloveshell.hh:58
typename Traits::BasisHandler BasisHandler
Definition: kirchhoffloveshell.hh:54
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:342
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:388
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:316
int order() const
Definition: kirchhoffloveshell.hh:154
A PreFE struct for Kirchhoff-Love shell elements.
Definition: kirchhoffloveshell.hh:31
double thickness
Definition: kirchhoffloveshell.hh:33
YoungsModulusAndPoissonsRatio material
Definition: kirchhoffloveshell.hh:32
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:434
double youngs_modulus
Definition: kirchhoffloveshell.hh:435
double nu
Definition: kirchhoffloveshell.hh:436
double thickness
Definition: kirchhoffloveshell.hh:437
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
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:32