version 0.4.1
nonlinearelastic.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 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 #include <dune/fufem/boundarypatch.hh>
14 #include <dune/geometry/quadraturerules.hh>
15 #include <dune/geometry/type.hh>
16 #include <dune/localfefunctions/cachedlocalBasis/cachedlocalBasis.hh>
17 #include <dune/localfefunctions/expressions/greenLagrangeStrains.hh>
18 #include <dune/localfefunctions/impl/standardLocalFunction.hh>
19 #include <dune/localfefunctions/manifolds/realTuple.hh>
20
31
32namespace Ikarus {
33
34template <typename PreFE, typename FE, typename PRE>
35class NonLinearElastic;
36
41template <Concepts::Material MAT>
43{
44 using Material = MAT;
46
47 template <typename PreFE, typename FE>
49};
50
60template <typename PreFE, typename FE, typename PRE>
61class NonLinearElastic : public ResultTypeBase<ResultTypes::PK2Stress, ResultTypes::PK2StressFull>
62{
63public:
65 using Basis = typename Traits::Basis;
66 using FlatBasis = typename Traits::FlatBasis;
68 using LocalView = typename Traits::LocalView;
69 using Geometry = typename Traits::Geometry;
70 using GridView = typename Traits::GridView;
71 using Element = typename Traits::Element;
72 using Material = PRE::Material;
73 using Pre = PRE;
74
75 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
76
77 template <typename ST>
78 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
79
80 static constexpr int myDim = Traits::mydim;
81 static constexpr int strainDim = myDim * (myDim + 1) / 2;
82 static constexpr auto strainType = StrainTags::greenLagrangian;
83 static constexpr auto stressType = StressTags::PK2;
84
85 template <template <typename, int, int> class RT>
87
88 template <typename ST = double>
89 using StrainType = Eigen::Vector<ST, strainDim>;
90
91 template <typename ST = double>
92 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
93
94 template <typename ST = double>
95 using KgType = Eigen::Matrix<ST, myDim, myDim>;
96
101 explicit NonLinearElastic(const Pre& pre)
102 : mat_{pre.material} {}
103
104protected:
108 void bindImpl() {
109 const auto& localView = underlying().localView();
110 const auto& element = localView.element();
111 auto& firstChild = localView.tree().child(0);
112 const auto& fe = firstChild.finiteElement();
113 geo_ = std::make_shared<const Geometry>(element.geometry());
114 numberOfNodes_ = fe.size();
115 order_ = 2 * (fe.localBasis().order());
116 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
117 if constexpr (requires { element.impl().getQuadratureRule(order_); })
118 if (element.impl().isTrimmed())
119 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
120 else
121 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
122 Dune::bindDerivatives(0, 1));
123 else
124 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
125 }
126
127public:
136 template <typename ScalarType = double>
137 auto displacementFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
138 const auto& d = par.globalSolution();
139 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
140 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
141 return uFunction;
142 }
143
152 template <typename ScalarType = double>
153 inline auto strainFunction(const Requirement& par, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
154 return Dune::greenLagrangeStrains(displacementFunction(par, dx));
155 }
156
165 template <typename ScalarType, int strainDim>
166 auto internalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
167 return material<ScalarType>().template storedEnergy<strainType>(strain);
168 }
169
179 template <typename ScalarType, int strainDim, bool voigt = true>
180 auto stress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
181 return material<ScalarType>().template stresses<strainType, voigt>(strain);
182 }
183
193 template <typename ScalarType, int strainDim, bool voigt = true>
194 auto materialTangent(const Eigen::Vector<ScalarType, strainDim>& strain) const {
195 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
196 }
197
198 const Geometry& geometry() const { return *geo_; }
199 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
200 [[nodiscard]] int order() const { return order_; }
201 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>& localBasis() const { return localBasis_; }
202
203 template <typename ScalarType = double>
204 decltype(auto) material() const {
206 return mat_.template rebind<ScalarType>();
207 else
208 return mat_;
209 }
210
211public:
217 template <template <typename, int, int> class RT>
218 requires(supportsResultType<RT>())
219 auto resultFunction() const {
220 return [&](const Eigen::Vector<double, strainDim>& strainInVoigt) {
221 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>) {
222 decltype(auto) mat = [&]() {
223 if constexpr (isSameResultType<RT, ResultTypes::PK2StressFull> and Material::isReduced)
224 return mat_.underlying();
225 else
226 return mat_;
227 }();
228 return RTWrapperType<RT>{mat.template stresses<strainType>(enlargeIfReduced<Material>(strainInVoigt))};
229 }
230 };
231 }
232
242 template <template <typename, int, int> class RT>
243 requires(supportsResultType<RT>())
245 Dune::PriorityTag<1>) const {
246 using namespace Dune::DerivativeDirections;
247
248 if constexpr (FE::isMixed())
249 return RTWrapperType<RT>{};
250 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>) {
251 const auto uFunction = displacementFunction(req);
252 const auto rFunction = resultFunction<RT>();
253 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
254 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
255
256 return rFunction(toVoigt(E));
257 }
258 }
259
260private:
261 //> CRTP
262 const auto& underlying() const { return static_cast<const FE&>(*this); }
263 auto& underlying() { return static_cast<FE&>(*this); }
264 std::shared_ptr<const Geometry> geo_;
265 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
266 Material mat_;
267 size_t numberOfNodes_{0};
268 int order_{};
269
270public:
281 template <typename ST>
282 auto geometricStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
283 const VectorXOptRef<ST>& dx = std::nullopt) const {
284 return [&](const FE::template KgType<ST>& kgIJ, const int I, const int J, const auto& gp) {
285 const auto geo = underlying().localView().element().geometry();
286 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
287 K.template block<FE::myDim, FE::myDim>(I * FE::myDim, J * FE::myDim) += kgIJ * intElement;
288 };
289 }
300 template <typename ST>
301 auto materialStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
302 const VectorXOptRef<ST>& dx = std::nullopt) const {
303 return [&](const StrainType<ST>& strain, const BopType<ST>& bopI, const BopType<ST>& bopJ, const int I, const int J,
304 const auto& gp) {
305 const auto C = materialTangent(strain);
306 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
307 K.template block<myDim, myDim>(I * myDim, J * myDim) += (bopI.transpose() * C * bopJ) * intElement;
308 };
309 }
310
322 template <typename ST>
323 auto internalForcesFunction(const Requirement& par, typename Traits::template VectorType<ST>& force,
324 const VectorXOptRef<ST>& dx = std::nullopt) const {
325 return [&](const StrainType<ST>& stresses, const BopType<ST>& bopI, const int I, const auto& gp) {
326 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
327 force.template segment<myDim>(myDim * I) += bopI.transpose() * stresses * intElement;
328 };
329 }
330
339 template <typename ST>
340 auto energyFunction(const Requirement& par, const VectorXOptRef<ST>& dx = std::nullopt) const {
341 return [&]() -> ST {
342 using namespace Dune::DerivativeDirections;
343 using namespace Dune;
344 ST energy = 0.0;
345 const auto eps = strainFunction(par, dx);
346 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
347 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
348 energy += internalEnergy(EVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
349 }
350 return energy;
351 };
352 }
353
354protected:
362 template <typename ScalarType>
363 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
364 typename Traits::template MatrixType<> K,
365 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
366 if constexpr (FE::isMixed())
367 return;
368 using namespace Dune::DerivativeDirections;
369 using namespace Dune;
370 const auto uFunction = displacementFunction(par, dx);
371 const auto eps = strainFunction(par, dx);
372 const auto kMFunction = materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
373 const auto kGFunction = geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
374 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
375 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
376 const auto stresses = stress(EVoigt);
377 for (size_t i = 0; i < numberOfNodes_; ++i) {
378 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
379 for (size_t j = 0; j < numberOfNodes_; ++j) {
380 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
381 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
382 kMFunction(EVoigt, bopI, bopJ, i, j, gp);
383 kGFunction(kgIJ, i, j, gp);
384 }
385 }
386 }
387 }
388
389 template <typename ScalarType>
391 const VectorXOptRef<ScalarType>& dx = std::nullopt) const -> ScalarType {
392 if constexpr (FE::isMixed())
393 return ScalarType{0.0};
394 return energyFunction(par, dx)();
395 }
396
397 template <typename ScalarType>
399 typename Traits::template VectorType<ScalarType> force,
400 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
401 if constexpr (FE::isMixed())
402 return;
403 using namespace Dune::DerivativeDirections;
404 using namespace Dune;
405 const auto eps = strainFunction(par, dx);
406 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
407
408 // Internal forces
409 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
410 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
411 const auto stresses = stress(EVoigt);
412 for (size_t i = 0; i < numberOfNodes_; ++i) {
413 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
414 fIntFunction(stresses, bopI, i, gp);
415 }
416 }
417 }
418};
419
426template <Concepts::Material MAT>
427auto nonLinearElastic(const MAT& mat) {
429
430 return pre;
431}
432
433} // namespace Ikarus
434
435#else
436 #error NonLinearElastic depends on dune-localfefunctions, which is not included
437#endif
Collection of fallback default functions.
Helper for the autodiff library.
Helper for transform between Dune linear algebra types and Eigen.
Definition of several material related enums.
Header file for types of loads in Ikarus finite element mechanics.
Header file for various EAS functions.
Definition of the LinearElastic class for finite element mechanics computations.
Material property functions and conversion utilities.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
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:64
auto nonLinearElastic(const MAT &mat)
A helper function to create a non-linear elastic pre finite element.
Definition: nonlinearelastic.hh:427
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
Definition: utils/dirichletvalues.hh:30
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: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
std::conditional_t< useFlat, FlatBasis, UntouchedBasis > Basis
Type of the basis version.
Definition: fetraits.hh:39
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
NonLinearElastic class represents a non-linear elastic finite element.
Definition: nonlinearelastic.hh:62
auto strainFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
The strain function for the given Requirement.
Definition: nonlinearelastic.hh:153
typename Traits::Basis Basis
Definition: nonlinearelastic.hh:65
typename Traits::Geometry Geometry
Definition: nonlinearelastic.hh:69
typename Traits::FlatBasis FlatBasis
Definition: nonlinearelastic.hh:66
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: nonlinearelastic.hh:92
auto energyFunction(const Requirement &par, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the internal energy at a given integration point and its index.
Definition: nonlinearelastic.hh:340
PRE::Material Material
Definition: nonlinearelastic.hh:72
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > & localBasis() const
Definition: nonlinearelastic.hh:201
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: nonlinearelastic.hh:78
typename Traits::GridView GridView
Definition: nonlinearelastic.hh:70
void bindImpl()
A helper function to bind the local view to the element.
Definition: nonlinearelastic.hh:108
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: nonlinearelastic.hh:398
const Geometry & geometry() const
Definition: nonlinearelastic.hh:198
static constexpr int myDim
Definition: nonlinearelastic.hh:80
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: nonlinearelastic.hh:390
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: nonlinearelastic.hh:244
static constexpr auto stressType
Definition: nonlinearelastic.hh:83
decltype(auto) material() const
Definition: nonlinearelastic.hh:204
NonLinearElastic(const Pre &pre)
Constructor for the NonLinearElastic class.
Definition: nonlinearelastic.hh:101
typename Traits::Element Element
Definition: nonlinearelastic.hh:71
PRE Pre
Definition: nonlinearelastic.hh:73
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: nonlinearelastic.hh:95
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: nonlinearelastic.hh:282
auto resultFunction() const
Get a lambda function that evaluates the requested result type for a given strain (in Voigt notation)...
Definition: nonlinearelastic.hh:219
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the material tangent for the given strain for the given Requirement.
Definition: nonlinearelastic.hh:194
static constexpr int strainDim
Definition: nonlinearelastic.hh:81
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 material part of the stiffness matrix (Ke + Ku) for a given ...
Definition: nonlinearelastic.hh:301
auto internalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: nonlinearelastic.hh:166
static constexpr auto strainType
Definition: nonlinearelastic.hh:82
Eigen::Vector< ST, strainDim > StrainType
Definition: nonlinearelastic.hh:89
typename Traits::LocalView LocalView
Definition: nonlinearelastic.hh:68
int order() const
Definition: nonlinearelastic.hh:200
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: nonlinearelastic.hh:75
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: nonlinearelastic.hh:323
auto stress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: nonlinearelastic.hh:180
auto displacementFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get the displacement function for the given Requirement.
Definition: nonlinearelastic.hh:137
size_t numberOfNodes() const
Definition: nonlinearelastic.hh:199
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Calculate the matrix associated with the given Requirement.
Definition: nonlinearelastic.hh:363
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:43
MAT Material
Definition: nonlinearelastic.hh:44
MAT material
Definition: nonlinearelastic.hh:45
static constexpr bool isMixed()
Definition: mixin.hh:110
Definition: utils/dirichletvalues.hh:32
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:622