version 0.4.8
linearelastic.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2026 The Ikarus Developers ikarus@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
10#pragma once
11
12#if HAVE_DUNE_LOCALFEFUNCTIONS
13
14 #include <optional>
15 #include <type_traits>
16
17 #include <dune/geometry/quadraturerules.hh>
18 #include <dune/localfefunctions/expressions/linearStrainsExpr.hh>
19 #include <dune/localfefunctions/impl/standardLocalFunction.hh>
20
29
30namespace Ikarus {
31
32template <typename PreFE, typename FE, typename PRE>
33class LinearElastic;
34
39template <Concepts::GeometricallyLinearMaterial MAT>
41{
42 using Material = MAT;
44
45 template <typename PreFE, typename FE>
47};
48
57template <typename PreFE, typename FE, typename PRE>
58class LinearElastic : public ResultTypeBase<ResultTypes::linearStress, ResultTypes::linearStressFull>
59{
60public:
63 using FlatBasis = typename Traits::FlatBasis;
65 using LocalView = typename Traits::LocalView;
66 using Geometry = typename Traits::Geometry;
67 using GridView = typename Traits::GridView;
68 using Element = typename Traits::Element;
69 using Material = PRE::Material;
70 using Pre = PRE;
71
72 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
73
74 template <typename ST>
75 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
76
77 static constexpr int myDim = Traits::mydim;
78 static constexpr int strainDim = myDim * (myDim + 1) / 2;
79 static constexpr auto strainType = StrainTags::linear;
80 static constexpr auto stressType = StressTags::linear;
81
82 template <template <typename, int, int> class RT>
84
85 template <typename ST = double>
86 using StrainType = Eigen::Vector<ST, strainDim>;
87
88 template <typename ST = double>
89 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
90
91 template <typename ST = double>
92 using KgType = Eigen::Matrix<ST, myDim, myDim>;
93
98 explicit LinearElastic(const Pre& pre)
99 : mat_{pre.material} {}
100
101protected:
105 void bindImpl() {
106 const auto& localView = underlying().localView();
107 const auto& element = localView.element();
108 auto& firstChild = localView.tree().child(0);
109 const auto& fe = firstChild.finiteElement();
110 geo_ = std::make_shared<const Geometry>(element.geometry());
111 numberOfNodes_ = fe.size();
112 order_ = fe.localBasis().order();
113 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
114 if (underlying().quadratureRule())
115 localBasis_.bind(underlying().quadratureRule().value(), Dune::bindDerivatives(0, 1));
116 else
117 localBasis_.bind(defaultQuadratureRule<myDim>(element, 2 * order_), Dune::bindDerivatives(0, 1));
118 }
119
120public:
129 template <typename ScalarType = double>
130 auto displacementFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
131 const auto& d = par.globalSolution();
132 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
133 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
134 return uFunction;
135 }
145 template <class ScalarType = double>
146 auto strainFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
147 return Dune::linearStrains(displacementFunction(par, dx));
148 }
149
158 template <typename ScalarType, int strainDim>
159 auto internalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
160 const auto C = materialTangent<ScalarType, strainDim>(strain);
161 return (0.5 * strain.dot(C * strain));
162 }
163
173 template <typename ScalarType, int strainDim, bool voigt = true>
174 auto stress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
175 const auto C = materialTangent<ScalarType, strainDim, voigt>(strain);
176 return (C * strain).eval();
177 }
178
188 template <typename ScalarType, int strainDim, bool voigt = true>
190 const Eigen::Vector<ScalarType, strainDim>& strain = Eigen::Vector<ScalarType, strainDim>::Zero()) const {
191 // Since that material is independent of the strains, a zero strain is passed here
192 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
193 }
194
195 const Geometry& geometry() const { return *geo_; }
196 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
197 [[nodiscard]] int order() const { return order_; }
198 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>& localBasis() const { return localBasis_; }
199
200 template <typename ScalarType = double>
201 decltype(auto) material() const {
203 return mat_.template rebind<ScalarType>();
204 else
205 return mat_;
206 }
207
208public:
220 template <typename M>
221 auto calculateStress(const M& mat, const Eigen::Vector<double, strainDim>& strainInVoigt) const {
222 return mat.template stresses<strainType>(enlargeIfReduced<Material>(strainInVoigt)).eval();
223 }
224
233 template <template <typename, int, int> class RT>
234 requires(LinearElastic::template supportsResultType<RT>())
236 Dune::PriorityTag<1>) const {
237 if constexpr (FE::isMixed())
238 return RTWrapperType<RT>{};
239
240 const auto eps = strainFunction(req);
241 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
242 decltype(auto) mat = [&]() {
243 if constexpr (isSameResultType<RT, ResultTypes::linearStressFull> and requires { mat_.underlying(); })
244 return mat_.underlying();
245 else
246 return mat_;
247 }();
248 const auto linearStress = calculateStress(mat, epsVoigt).eval();
249 RTWrapperType<RT> resultWrapper{};
250 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or
251 isSameResultType<RT, ResultTypes::linearStressFull>)
252 resultWrapper = linearStress;
253 return resultWrapper;
254 }
255
256private:
257 //> CRTP
258 const auto& underlying() const { return static_cast<const FE&>(*this); }
259 auto& underlying() { return static_cast<FE&>(*this); }
260
261 std::shared_ptr<const Geometry> geo_;
262 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
263 Material mat_;
264 size_t numberOfNodes_{0};
265 int order_{};
266
267public:
279 template <typename ST>
280 auto geometricStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
281 const VectorXOptRef<ST>& dx = std::nullopt) const {
282 return [&](const FE::template KgType<ST>& kgIJ, const int I, const int J, const auto& gp) {
283 const auto geo = underlying().localView().element().geometry();
284 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
285 K.template block<FE::myDim, FE::myDim>(I * FE::myDim, J * FE::myDim) += kgIJ * intElement;
286 };
287 }
298 template <typename ST>
299 auto materialStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
300 const VectorXOptRef<ST>& dx = std::nullopt) const {
301 return [&](const StrainType<ST>& strain, const BopType<ST>& bopI, const BopType<ST>& bopJ, const int I, const int J,
302 const auto& gp) {
303 const auto C = materialTangent(strain);
304 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
305 K.template block<myDim, myDim>(I * myDim, J * myDim) += (bopI.transpose() * C * bopJ) * intElement;
306 };
307 }
308
320 template <typename ST>
321 auto internalForcesFunction(const Requirement& par, typename Traits::template VectorType<ST>& force,
322 const VectorXOptRef<ST>& dx = std::nullopt) const {
323 return [&](const StrainType<ST>& stresses, const BopType<ST>& bopI, const int I, const auto& gp) {
324 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
325 force.template segment<myDim>(myDim * I) += bopI.transpose() * stresses * intElement;
326 };
327 }
328
337 template <typename ScalarType>
338 auto energyFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
339 return [&]() -> ScalarType {
340 using namespace Dune::DerivativeDirections;
341 using namespace Dune;
342 ScalarType energy = 0.0;
343 const auto eps = strainFunction(par, dx);
344 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
345 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
346 energy += internalEnergy(epsVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
347 }
348 return energy;
349 };
350 }
351
352protected:
353 template <typename ScalarType>
354 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
355 typename Traits::template MatrixType<> K,
356 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
357 if constexpr (FE::isMixed()) {
358 return;
359 }
360 using namespace Dune::DerivativeDirections;
361 using namespace Dune;
362 const auto eps = strainFunction(par, dx);
363 const auto kFunction = materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
364
365 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
366 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
367 for (size_t i = 0; i < numberOfNodes_; ++i) {
368 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
369 for (size_t j = 0; j < numberOfNodes_; ++j) {
370 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
371 kFunction(epsVoigt, bopI, bopJ, i, j, gp);
372 }
373 }
374 }
375 }
376
377 template <typename ScalarType>
379 const VectorXOptRef<ScalarType>& dx = std::nullopt) const -> ScalarType {
380 if constexpr (FE::isMixed())
381 return ScalarType{0.0};
382 return energyFunction(par, dx)();
383 }
384
385 template <typename ScalarType>
387 typename Traits::template VectorType<ScalarType> force,
388 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
389 if constexpr (FE::isMixed())
390 return;
391 using namespace Dune::DerivativeDirections;
392 using namespace Dune;
393 const auto eps = strainFunction(par, dx);
394 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
395
396 // Internal forces
397 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
398 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
399 const auto stresses = stress(epsVoigt);
400 for (size_t i = 0; i < numberOfNodes_; ++i) {
401 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
402 fIntFunction(stresses, bopI, i, gp);
403 }
404 }
405 }
406};
407
414template <Concepts::GeometricallyLinearMaterial MAT>
415auto linearElastic(const MAT& mat) {
416 LinearElasticPre<MAT> pre(mat);
417
418 return pre;
419}
420} // namespace Ikarus
421
422#else
423 #error LinearElastic depends on dune-localfefunctions, which is not included
424#endif
Free function for creating a tensor product quadrature rule and other Dune::Quadrature rule utils.
Definitions of ResultTypes used for finite element results.
Material property functions and conversion utilities.
Definition of several material related enums.
Header file for various EAS functions.
Definition of the LinearElastic class for finite element mechanics computations.
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
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:39
auto linearElastic(const MAT &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:415
Definition: utils/dirichletvalues.hh:36
FE class is a base class for all finite elements.
Definition: febase.hh:79
static constexpr int myDim
Definition: febase.hh:95
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
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:165
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:278
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
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
LinearElastic class represents a linear elastic finite element.
Definition: linearelastic.hh:59
void bindImpl()
A helper function to bind the local view to the element.
Definition: linearelastic.hh:105
auto internalForcesFunction(const Requirement &par, typename Traits::template VectorType< ST > &force, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the internal force vector for a given strain,...
Definition: linearelastic.hh:321
auto energyFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get a lambda function that evaluates the internal energy at a given integration point and its index.
Definition: linearelastic.hh:338
static constexpr int myDim
Definition: linearelastic.hh:77
Eigen::Vector< ST, strainDim > StrainType
Definition: linearelastic.hh:86
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: linearelastic.hh:89
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: linearelastic.hh:354
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: linearelastic.hh:92
typename Traits::Element Element
Definition: linearelastic.hh:68
static constexpr auto strainType
Definition: linearelastic.hh:79
PRE Pre
Definition: linearelastic.hh:70
static constexpr auto stressType
Definition: linearelastic.hh:80
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: linearelastic.hh:378
size_t numberOfNodes() const
Definition: linearelastic.hh:196
int order() const
Definition: linearelastic.hh:197
auto materialStiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the elastic stiffness matrix for a given strain,...
Definition: linearelastic.hh:299
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain=Eigen::Vector< ScalarType, strainDim >::Zero()) const
Get the linear elastic material tangent for the given strain for the given Requirement.
Definition: linearelastic.hh:189
auto calculateAtImpl(const Requirement &req, const Dune::FieldVector< double, Traits::mydim > &local, Dune::PriorityTag< 1 >) const
Calculates a requested result at a specific local position.
Definition: linearelastic.hh:235
auto geometricStiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the geometric part of the stiffness matrix (Kg) for a given inte...
Definition: linearelastic.hh:280
typename Traits::GridView GridView
Definition: linearelastic.hh:67
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:62
typename Traits::LocalView LocalView
Definition: linearelastic.hh:65
auto calculateStress(const M &mat, const Eigen::Vector< double, strainDim > &strainInVoigt) const
Get the linear stress for a given strain (in Voigt notation).
Definition: linearelastic.hh:221
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: linearelastic.hh:75
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:63
auto displacementFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Gets the displacement function for the given Requirement and optional displacement vector.
Definition: linearelastic.hh:130
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:98
PRE::Material Material
Definition: linearelastic.hh:69
static constexpr int strainDim
Definition: linearelastic.hh:78
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:72
decltype(auto) material() const
Definition: linearelastic.hh:201
typename Traits::Geometry Geometry
Definition: linearelastic.hh:66
auto internalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: linearelastic.hh:159
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > & localBasis() const
Definition: linearelastic.hh:198
auto strainFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Gets the strain function for the given Requirement and optional di splacement vector.
Definition: linearelastic.hh:146
const Geometry & geometry() const
Definition: linearelastic.hh:195
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: linearelastic.hh:386
auto stress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: linearelastic.hh:174
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:41
MAT Material
Definition: linearelastic.hh:42
MAT material
Definition: linearelastic.hh:43
static constexpr bool isMixed()
Definition: mixin.hh:110
Definition: utils/dirichletvalues.hh:38
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:625
Header file for material models in Ikarus finite element mechanics.