version 0.4.1
linearelastic.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers mueller@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
28
29namespace Ikarus {
30
31template <typename PreFE, typename FE, typename PRE>
32class LinearElastic;
33
38template <Concepts::GeometricallyLinearMaterial MAT>
40{
41 using Material = MAT;
43
44 template <typename PreFE, typename FE>
46};
47
56template <typename PreFE, typename FE, typename PRE>
57class LinearElastic : public ResultTypeBase<ResultTypes::linearStress, ResultTypes::linearStressFull>
58{
59public:
62 using FlatBasis = typename Traits::FlatBasis;
64 using LocalView = typename Traits::LocalView;
65 using Geometry = typename Traits::Geometry;
66 using GridView = typename Traits::GridView;
67 using Element = typename Traits::Element;
68 using Material = PRE::Material;
69 using Pre = PRE;
70
71 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
72
73 template <typename ST>
74 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
75
76 static constexpr int myDim = Traits::mydim;
77 static constexpr int strainDim = myDim * (myDim + 1) / 2;
78 static constexpr auto strainType = StrainTags::linear;
79 static constexpr auto stressType = StressTags::linear;
80 static constexpr bool hasEAS = FE::template hasEAS<EAS::LinearStrain>;
81
82 template <template <typename, int, int> class RT>
84
85 template <typename ST>
86 using StrainType = Eigen::Vector<ST, strainDim>;
87
88 template <typename ST>
89 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
90
91 template <typename ST>
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_ = 2 * (fe.localBasis().order());
113 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
114 if constexpr (requires { element.impl().getQuadratureRule(order_); })
115 if (element.impl().isTrimmed())
116 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
117 else
118 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
119 Dune::bindDerivatives(0, 1));
120 else
121 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
122 }
123
124public:
133 template <typename ScalarType = double>
134 auto displacementFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
135 const auto& d = par.globalSolution();
136 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
137 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
138 return uFunction;
139 }
149 template <class ScalarType = double>
150 auto strainFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
151 return Dune::linearStrains(displacementFunction(par, dx));
152 }
153
162 template <typename ScalarType, int strainDim>
163 auto internalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
164 const auto C = materialTangent<ScalarType, strainDim>(strain);
165 return (0.5 * strain.dot(C * strain));
166 }
167
177 template <typename ScalarType, int strainDim, bool voigt = true>
178 auto stress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
179 const auto C = materialTangent<ScalarType, strainDim, voigt>(strain);
180 return (C * strain).eval();
181 }
182
192 template <typename ScalarType, int strainDim, bool voigt = true>
194 const Eigen::Vector<ScalarType, strainDim>& strain = Eigen::Vector<ScalarType, strainDim>::Zero()) const {
195 // Since that material is independent of the strains, a zero strain is passed here
196 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
197 }
198
199 const Geometry& geometry() const { return *geo_; }
200 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
201 [[nodiscard]] int order() const { return order_; }
202 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>& localBasis() const { return localBasis_; }
203
204 template <typename ScalarType = double>
205 decltype(auto) material() const {
207 return mat_.template rebind<ScalarType>();
208 else
209 return mat_;
210 }
211
212public:
218 template <template <typename, int, int> class RT>
219 requires(LinearElastic::template supportsResultType<RT>())
220 auto resultFunction() const {
221 return [&](const Eigen::Vector<double, strainDim>& strainInVoigt) {
222 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or
223 isSameResultType<RT, ResultTypes::linearStressFull>) {
224 decltype(auto) mat = [&]() {
225 if constexpr (isSameResultType<RT, ResultTypes::linearStressFull> and requires { mat_.underlying(); })
226 return mat_.underlying();
227 else
228 return mat_;
229 }();
230 return RTWrapperType<RT>{mat.template stresses<StrainTags::linear>(enlargeIfReduced<Material>(strainInVoigt))};
231 }
232 };
233 }
234
243 template <template <typename, int, int> class RT>
244 requires(LinearElastic::template supportsResultType<RT>())
246 Dune::PriorityTag<1>) const {
247 if constexpr (hasEAS)
248 return RTWrapperType<RT>{};
249 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or
250 isSameResultType<RT, ResultTypes::linearStressFull>) {
251 const auto rFunction = resultFunction<RT>();
252 const auto eps = strainFunction(req);
253 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
254 return rFunction(epsVoigt);
255 }
256 }
257
258private:
259 //> CRTP
260 const auto& underlying() const { return static_cast<const FE&>(*this); }
261 auto& underlying() { return static_cast<FE&>(*this); }
262
263 std::shared_ptr<const Geometry> geo_;
264 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
265 Material mat_;
266 size_t numberOfNodes_{0};
267 int order_{};
268
269public:
280 template <typename ST>
281 auto stiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
282 const VectorXOptRef<ST>& dx = std::nullopt) const {
283 return [&](const StrainType<ST>& strain, const BopType<ST>& bopI, const BopType<ST>& bopJ, const KgType<ST>& kgIJ,
284 const int I, const int J, const auto& gp) {
285 const auto C = materialTangent(strain);
286 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
287 K.template block<myDim, myDim>(I * myDim, J * myDim) += (bopI.transpose() * C * bopJ) * intElement;
288 };
289 }
290
302 template <typename ST>
303 auto internalForcesFunction(const Requirement& par, typename Traits::template VectorType<ST>& force,
304 const VectorXOptRef<ST>& dx = std::nullopt) const {
305 return [&](const StrainType<ST>& strain, const BopType<ST>& bopI, const int I, const auto& gp) {
306 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
307 auto stresses = stress(strain);
308 force.template segment<myDim>(myDim * I) += bopI.transpose() * stresses * intElement;
309 };
310 }
311
320 template <typename ScalarType>
321 auto energyFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
322 return [&]() -> ScalarType {
323 using namespace Dune::DerivativeDirections;
324 using namespace Dune;
325 ScalarType energy = 0.0;
326 const auto eps = strainFunction(par, dx);
327 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
328 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
329 energy += internalEnergy(epsVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
330 }
331 return energy;
332 };
333 }
334
335protected:
336 template <typename ScalarType>
337 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
338 typename Traits::template MatrixType<> K,
339 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
340 if constexpr (hasEAS)
341 return;
342 using namespace Dune::DerivativeDirections;
343 using namespace Dune;
344 const auto eps = strainFunction(par, dx);
345 const auto kFunction = stiffnessMatrixFunction<ScalarType>(par, K, dx);
346 const auto& kgIJ = KgType<ScalarType>::Zero().eval();
347
348 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
349 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
350 for (size_t i = 0; i < numberOfNodes_; ++i) {
351 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
352 for (size_t j = 0; j < numberOfNodes_; ++j) {
353 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
354 kFunction(epsVoigt, bopI, bopJ, kgIJ, i, j, gp);
355 }
356 }
357 }
358 }
359
360 template <typename ScalarType>
362 const VectorXOptRef<ScalarType>& dx = std::nullopt) const -> ScalarType {
363 if constexpr (hasEAS)
364 return ScalarType{0.0};
365 return energyFunction(par, dx)();
366 }
367
368 template <typename ScalarType>
370 typename Traits::template VectorType<ScalarType> force,
371 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
372 if constexpr (hasEAS)
373 return;
374 using namespace Dune::DerivativeDirections;
375 using namespace Dune;
376 const auto eps = strainFunction(par, dx);
377 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
378
379 // Internal forces
380 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
381 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
382 for (size_t i = 0; i < numberOfNodes_; ++i) {
383 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
384 fIntFunction(epsVoigt, bopI, i, gp);
385 }
386 }
387 }
388};
389
396template <typename MAT>
397auto linearElastic(const MAT& mat) {
398 LinearElasticPre<MAT> pre(mat);
399
400 return pre;
401}
402} // namespace Ikarus
403
404#else
405 #error LinearElastic depends on dune-localfefunctions, which is not included
406#endif
Definition of the LinearElastic class for finite element mechanics computations.
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: assemblermanipulatorbuildingblocks.hh:22
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:64
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:49
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:38
auto linearElastic(const MAT &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:397
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:223
const SolutionVectorType & globalSolution() const
Get the global solution vector.
Definition: ferequirements.hh:313
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
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:58
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:303
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:321
static constexpr int myDim
Definition: linearelastic.hh:76
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:337
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: linearelastic.hh:92
typename Traits::Element Element
Definition: linearelastic.hh:67
auto resultFunction() const
Get a lambda function that evaluates the requested result type for a given strain (in Voigt notation)...
Definition: linearelastic.hh:220
static constexpr auto strainType
Definition: linearelastic.hh:78
auto stiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the stiffness matrix for a given strain, integration point and i...
Definition: linearelastic.hh:281
PRE Pre
Definition: linearelastic.hh:69
static constexpr auto stressType
Definition: linearelastic.hh:79
static constexpr bool hasEAS
Definition: linearelastic.hh:80
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: linearelastic.hh:361
size_t numberOfNodes() const
Definition: linearelastic.hh:200
int order() const
Definition: linearelastic.hh:201
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:193
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:245
typename Traits::GridView GridView
Definition: linearelastic.hh:66
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:61
typename Traits::LocalView LocalView
Definition: linearelastic.hh:64
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: linearelastic.hh:74
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:62
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:134
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:98
PRE::Material Material
Definition: linearelastic.hh:68
static constexpr int strainDim
Definition: linearelastic.hh:77
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:71
decltype(auto) material() const
Definition: linearelastic.hh:205
typename Traits::Geometry Geometry
Definition: linearelastic.hh:65
auto internalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: linearelastic.hh:163
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > & localBasis() const
Definition: linearelastic.hh:202
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:150
const Geometry & geometry() const
Definition: linearelastic.hh:199
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: linearelastic.hh:369
auto stress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: linearelastic.hh:178
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:40
MAT Material
Definition: linearelastic.hh:41
MAT material
Definition: linearelastic.hh:42
Definition: utils/dirichletvalues.hh:32
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:615
Header file for material models in Ikarus finite element mechanics.