version 0.4.8
nonlinearelastic.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 #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
34
35namespace Ikarus {
36
37template <typename PreFE, typename FE, typename PRE>
38class NonLinearElastic;
39
44template <Concepts::Material MAT>
46{
47 using Material = MAT;
49
50 template <typename PreFE, typename FE>
52};
53
63template <typename PreFE, typename FE, typename PRE>
64class NonLinearElastic : public ResultTypeBase<ResultTypes::PK2Stress, ResultTypes::PK2StressFull,
65 ResultTypes::kirchhoffStress, ResultTypes::cauchyStress>
66{
67public:
69 using Basis = typename Traits::Basis;
70 using FlatBasis = typename Traits::FlatBasis;
72 using LocalView = typename Traits::LocalView;
73 using Geometry = typename Traits::Geometry;
74 using GridView = typename Traits::GridView;
75 using Element = typename Traits::Element;
76 using Material = PRE::Material;
77 using Pre = PRE;
78
79 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
80
81 template <typename ST>
82 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
83
84 static constexpr int myDim = Traits::mydim;
85 static constexpr int strainDim = myDim * (myDim + 1) / 2;
86 static constexpr auto strainType = StrainTags::greenLagrangian;
87 static constexpr auto stressType = StressTags::PK2;
88
89 template <template <typename, int, int> class RT>
91
92 template <typename ST = double>
93 using StrainType = Eigen::Vector<ST, strainDim>;
94
95 template <typename ST = double>
96 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
97
98 template <typename ST = double>
99 using KgType = Eigen::Matrix<ST, myDim, myDim>;
100
105 explicit NonLinearElastic(const Pre& pre)
106 : mat_{pre.material} {}
107
108protected:
112 void bindImpl() {
113 const auto& localView = underlying().localView();
114 const auto& element = localView.element();
115 auto& firstChild = localView.tree().child(0);
116 const auto& fe = firstChild.finiteElement();
117 geo_ = std::make_shared<const Geometry>(element.geometry());
118 numberOfNodes_ = fe.size();
119 order_ = fe.localBasis().order();
120 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
121 if (underlying().quadratureRule())
122 localBasis_.bind(underlying().quadratureRule().value(), Dune::bindDerivatives(0, 1));
123 else
124 localBasis_.bind(defaultQuadratureRule<myDim>(element, 2 * 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:
223 template <typename M>
224 auto calculateStress(const M& mat, const Eigen::Vector<double, strainDim>& strainInVoigt) const {
225 return mat.template stresses<strainType>(enlargeIfReduced<Material>(strainInVoigt)).eval();
226 }
227
237 template <template <typename, int, int> class RT>
238 requires(supportsResultType<RT>())
240 Dune::PriorityTag<1>) const {
241 if constexpr (FE::isMixed())
242 return RTWrapperType<RT>{};
243
244 using namespace Dune::DerivativeDirections;
247 const auto uFunction = displacementFunction(req);
248 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
249 const auto F = transformStrain<StrainTags::displacementGradient, StrainTags::deformationGradient>(H).eval();
250 const auto E = transformStrain<StrainTags::displacementGradient, StrainTags::greenLagrangian>(H).eval();
251 decltype(auto) mat = [&]() {
252 if constexpr (isSameResultType<RT, ResultTypes::PK2StressFull> and Material::isReduced)
253 return mat_.underlying();
254 else
255 return mat_;
256 }();
257
258 const auto PK2Stress = calculateStress(mat, toVoigt(E)).eval();
259 RTWrapperType<RT> resultWrapper{};
260
261 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>)
262 resultWrapper = PK2Stress;
263 else if constexpr (isSameResultType<RT, ResultTypes::kirchhoffStress> or
264 isSameResultType<RT, ResultTypes::cauchyStress>) {
265 const auto PK2StressMat = fromVoigt(PK2Stress, false).eval();
266 constexpr StressTags toStress =
267 isSameResultType<RT, ResultTypes::kirchhoffStress> ? StressTags::Kirchhoff : StressTags::Cauchy;
268 resultWrapper = toVoigt(transformStress<StressTags::PK2, toStress>(PK2StressMat, F), false);
269 }
270 return resultWrapper;
271 }
272
273private:
274 //> CRTP
275 const auto& underlying() const { return static_cast<const FE&>(*this); }
276 auto& underlying() { return static_cast<FE&>(*this); }
277 std::shared_ptr<const Geometry> geo_;
278 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
279 Material mat_;
280 size_t numberOfNodes_{0};
281 int order_{};
282
283public:
294 template <typename ST>
295 auto geometricStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
296 const VectorXOptRef<ST>& dx = std::nullopt) const {
297 return [&](const FE::template KgType<ST>& kgIJ, const int I, const int J, const auto& gp) {
298 const auto geo = underlying().localView().element().geometry();
299 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
300 K.template block<FE::myDim, FE::myDim>(I * FE::myDim, J * FE::myDim) += kgIJ * intElement;
301 };
302 }
313 template <typename ST>
314 auto materialStiffnessMatrixFunction(const Requirement& par, typename Traits::template MatrixType<ST>& K,
315 const VectorXOptRef<ST>& dx = std::nullopt) const {
316 return [&](const StrainType<ST>& strain, const BopType<ST>& bopI, const BopType<ST>& bopJ, const int I, const int J,
317 const auto& gp) {
318 const auto C = materialTangent(strain);
319 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
320 K.template block<myDim, myDim>(I * myDim, J * myDim) += (bopI.transpose() * C * bopJ) * intElement;
321 };
322 }
323
335 template <typename ST>
336 auto internalForcesFunction(const Requirement& par, typename Traits::template VectorType<ST>& force,
337 const VectorXOptRef<ST>& dx = std::nullopt) const {
338 return [&](const StrainType<ST>& stresses, const BopType<ST>& bopI, const int I, const auto& gp) {
339 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
340 force.template segment<myDim>(myDim * I) += bopI.transpose() * stresses * intElement;
341 };
342 }
343
352 template <typename ST>
353 auto energyFunction(const Requirement& par, const VectorXOptRef<ST>& dx = std::nullopt) const {
354 return [&]() -> ST {
355 using namespace Dune::DerivativeDirections;
356 using namespace Dune;
357 ST energy = 0.0;
358 const auto eps = strainFunction(par, dx);
359 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
360 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
361 energy += internalEnergy(EVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
362 }
363 return energy;
364 };
365 }
366
367protected:
375 template <typename ScalarType>
376 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
377 typename Traits::template MatrixType<> K,
378 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
379 if constexpr (FE::isMixed())
380 return;
381 using namespace Dune::DerivativeDirections;
382 using namespace Dune;
383 const auto uFunction = displacementFunction(par, dx);
384 const auto eps = strainFunction(par, dx);
385 const auto kMFunction = materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
386 const auto kGFunction = geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
387 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
388 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
389 const auto stresses = stress(EVoigt);
390 for (size_t i = 0; i < numberOfNodes_; ++i) {
391 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
392 for (size_t j = 0; j < numberOfNodes_; ++j) {
393 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
394 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
395 kMFunction(EVoigt, bopI, bopJ, i, j, gp);
396 kGFunction(kgIJ, i, j, gp);
397 }
398 }
399 }
400 }
401
402 template <typename ScalarType>
404 const VectorXOptRef<ScalarType>& dx = std::nullopt) const -> ScalarType {
405 if constexpr (FE::isMixed())
406 return ScalarType{0.0};
407 return energyFunction(par, dx)();
408 }
409
410 template <typename ScalarType>
412 typename Traits::template VectorType<ScalarType> force,
413 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
414 if constexpr (FE::isMixed())
415 return;
416 using namespace Dune::DerivativeDirections;
417 using namespace Dune;
418 const auto eps = strainFunction(par, dx);
419 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
420
421 // Internal forces
422 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
423 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
424 const auto stresses = stress(EVoigt);
425 for (size_t i = 0; i < numberOfNodes_; ++i) {
426 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
427 fIntFunction(stresses, bopI, i, gp);
428 }
429 }
430 }
431};
432
439template <Concepts::Material MAT>
440auto nonLinearElastic(const MAT& mat) {
442
443 return pre;
444}
445
446} // namespace Ikarus
447
448#else
449 #error NonLinearElastic depends on dune-localfefunctions, which is not included
450#endif
Helper for transform between Dune linear algebra types and Eigen.
Helper for the autodiff library.
Free function for creating a tensor product quadrature rule and other Dune::Quadrature rule utils.
Collection of fallback default functions.
Material property functions and conversion utilities.
Header file for types of loads in Ikarus finite element mechanics.
Implementation of transformations for different strain measures.
Implementation of transformations for different stress measures.
Definition of several material related enums.
Header file for various EAS functions.
Definition of the LinearElastic class for finite element mechanics computations.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
StressTags
A strongly typed enum class representing the type of the computed stresses.
Definition: tags.hh:26
auto fromVoigt(const Eigen::Matrix< ST, size, 1, Options, maxSize, 1 > &EVoigt, bool isStrain=true)
Converts a vector given in Voigt notation to a matrix.
Definition: tensorutils.hh:296
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
auto transformStress(const Eigen::MatrixBase< DerivedS > &sRaw, const Eigen::MatrixBase< DerivedF > &F)
Transform stress measures from one type to another.
Definition: stressconversions.hh:142
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:65
auto nonLinearElastic(const MAT &mat)
A helper function to create a non-linear elastic pre finite element.
Definition: nonlinearelastic.hh:440
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:50
auto transformStrain(const Eigen::MatrixBase< Derived > &eRaw)
Transform strain from one type to another.
Definition: strainconversions.hh:119
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:39
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
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:66
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:69
typename Traits::Geometry Geometry
Definition: nonlinearelastic.hh:73
typename Traits::FlatBasis FlatBasis
Definition: nonlinearelastic.hh:70
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: nonlinearelastic.hh:96
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:353
PRE::Material Material
Definition: nonlinearelastic.hh:76
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:82
typename Traits::GridView GridView
Definition: nonlinearelastic.hh:74
void bindImpl()
A helper function to bind the local view to the element.
Definition: nonlinearelastic.hh:112
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: nonlinearelastic.hh:411
const Geometry & geometry() const
Definition: nonlinearelastic.hh:198
static constexpr int myDim
Definition: nonlinearelastic.hh:84
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: nonlinearelastic.hh:403
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:239
static constexpr auto stressType
Definition: nonlinearelastic.hh:87
decltype(auto) material() const
Definition: nonlinearelastic.hh:204
NonLinearElastic(const Pre &pre)
Constructor for the NonLinearElastic class.
Definition: nonlinearelastic.hh:105
typename Traits::Element Element
Definition: nonlinearelastic.hh:75
PRE Pre
Definition: nonlinearelastic.hh:77
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: nonlinearelastic.hh:99
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:295
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:85
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:314
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:86
auto calculateStress(const M &mat, const Eigen::Vector< double, strainDim > &strainInVoigt) const
Get the PK2 stress for a given strain (in Voigt notation).
Definition: nonlinearelastic.hh:224
Eigen::Vector< ST, strainDim > StrainType
Definition: nonlinearelastic.hh:93
typename Traits::LocalView LocalView
Definition: nonlinearelastic.hh:72
int order() const
Definition: nonlinearelastic.hh:200
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: nonlinearelastic.hh:79
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:336
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:376
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:46
MAT Material
Definition: nonlinearelastic.hh:47
MAT material
Definition: nonlinearelastic.hh:48
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