version 0.4.2
nonlinearelastic.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
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
30
31namespace Ikarus {
32
33template <typename PreFE, typename FE, typename PRE>
34class NonLinearElastic;
35
40template <typename MAT>
42{
43 using Material = MAT;
45
46 template <typename PreFE, typename FE>
48};
49
59template <typename PreFE, typename FE, typename PRE>
61{
62public:
64 using Basis = typename Traits::Basis;
65 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 static constexpr int myDim = Traits::mydim;
78 static constexpr auto strainType = StrainTags::greenLagrangian;
79
84 explicit NonLinearElastic(const Pre& pre)
85 : mat_{pre.material} {}
86
87protected:
91 void bindImpl() {
92 const auto& localView = underlying().localView();
93 const auto& element = localView.element();
94 auto& firstChild = localView.tree().child(0);
95 const auto& fe = firstChild.finiteElement();
96 geo_ = std::make_shared<const Geometry>(element.geometry());
97 numberOfNodes_ = fe.size();
98 order_ = 2 * (fe.localBasis().order());
99 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
100 if constexpr (requires { element.impl().getQuadratureRule(order_); })
101 if (element.impl().isTrimmed())
102 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
103 else
104 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
105 Dune::bindDerivatives(0, 1));
106 else
107 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
108 }
109
110public:
119 template <typename ScalarType = double>
121 const Requirement& par,
122 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
123 const auto& d = par.globalSolution();
124 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
125 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
126 return uFunction;
127 }
128
137 template <typename ScalarType = double>
138 inline auto strainFunction(
139 const Requirement& par,
140 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
141 return Dune::greenLagrangeStrains(displacementFunction(par, dx));
142 }
143
153 template <typename ScalarType, int strainDim, bool voigt = true>
154 auto materialTangent(const Eigen::Vector<ScalarType, strainDim>& strain) const {
155 if constexpr (std::is_same_v<ScalarType, double>)
156 return mat_.template tangentModuli<strainType, voigt>(strain);
157 else {
158 decltype(auto) matAD = mat_.template rebind<ScalarType>();
159 return matAD.template tangentModuli<strainType, voigt>(strain);
160 }
161 }
162
171 template <typename ScalarType, int strainDim>
172 auto getInternalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
173 if constexpr (std::is_same_v<ScalarType, double>)
174 return mat_.template storedEnergy<strainType>(strain);
175 else {
176 decltype(auto) matAD = mat_.template rebind<ScalarType>();
177 return matAD.template storedEnergy<strainType>(strain);
178 }
179 }
180
190 template <typename ScalarType, int strainDim, bool voigt = true>
191 auto getStress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
192 if constexpr (std::is_same_v<ScalarType, double>)
193 return mat_.template stresses<strainType, voigt>(strain);
194 else {
195 decltype(auto) matAD = mat_.template rebind<ScalarType>();
196 return matAD.template stresses<strainType, voigt>(strain);
197 }
198 }
199
200 const Geometry& geometry() const { return *geo_; }
201 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
202 [[nodiscard]] int order() const { return order_; }
203
204 using SupportedResultTypes = std::tuple<decltype(makeRT<ResultTypes::PK2Stress>())>;
205
206private:
207 template <template <typename, int, int> class RT>
208 static consteval bool canProvideResultType() {
209 return isSameResultType<RT, ResultTypes::PK2Stress>;
210 }
211
212public:
222 template <template <typename, int, int> class RT>
223 requires(canProvideResultType<RT>())
225 Dune::PriorityTag<1>) const {
226 using namespace Dune::DerivativeDirections;
227
229 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress>) {
230 const auto uFunction = displacementFunction(req);
231 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
232 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
233
234 return RTWrapper{mat_.template stresses<StrainTags::greenLagrangian>(toVoigt(E))};
235 }
236 }
237
238private:
239 //> CRTP
240 const auto& underlying() const { return static_cast<const FE&>(*this); }
241 auto& underlying() { return static_cast<FE&>(*this); }
242 std::shared_ptr<const Geometry> geo_;
243 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
244 Material mat_;
245 size_t numberOfNodes_{0};
246 int order_{};
247
248protected:
256 template <typename ScalarType>
258 const Requirement& par, const MatrixAffordance& affordance, typename Traits::template MatrixType<> K,
259 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
260 using namespace Dune::DerivativeDirections;
261 using namespace Dune;
262 const auto uFunction = displacementFunction(par, dx);
263 const auto eps = strainFunction(par, dx);
264 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
265 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
266 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
267 const auto u = (uFunction.evaluate(gpIndex, on(gridElement))).eval();
268 const auto C = materialTangent(EVoigt);
269
270 const auto stresses = getStress(EVoigt);
271 for (size_t i = 0; i < numberOfNodes_; ++i) {
272 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
273 for (size_t j = 0; j < numberOfNodes_; ++j) {
274 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
275 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
276 K.template block<myDim, myDim>(i * myDim, j * myDim) += (bopI.transpose() * C * bopJ + kgIJ) * intElement;
277 }
278 }
279 }
280 }
281
282 template <typename ScalarType>
284 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx =
285 std::nullopt) const -> ScalarType {
286 using namespace Dune::DerivativeDirections;
287 using namespace Dune;
288
289 const auto eps = strainFunction(par, dx);
290 const auto& lambda = par.parameter();
291 ScalarType energy = 0.0;
292
293 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
294 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
295 const auto internalEnergy = getInternalEnergy(EVoigt);
296 energy += internalEnergy * geo_->integrationElement(gp.position()) * gp.weight();
297 }
298
299 return energy;
300 }
301
302 template <typename ScalarType>
304 const Requirement& par, VectorAffordance affordance, typename Traits::template VectorType<ScalarType> force,
305 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
306 using namespace Dune::DerivativeDirections;
307 using namespace Dune;
308 const auto eps = strainFunction(par, dx);
309
310 // Internal forces
311 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
312 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
313 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
314 const auto stresses = getStress(EVoigt);
315 for (size_t i = 0; i < numberOfNodes_; ++i) {
316 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
317 force.template segment<myDim>(myDim * i) += bopI.transpose() * stresses * intElement;
318 }
319 }
320 }
321};
322
329template <typename MAT>
330auto nonLinearElastic(const MAT& mat) {
332
333 return pre;
334}
335
336} // namespace Ikarus
337
338#else
339 #error NonLinearElastic depends on dune-localfefunctions, which is not included
340#endif
Helper for the autodiff library.
Collection of fallback default functions.
Helper for transform between Dune linear algebra types and Eigen.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
Definition of the LinearElastic class for finite element mechanics computations.
Definition of several material related enums.
Header file for types of loads in Ikarus finite element mechanics.
Material property functions and conversion utilities.
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:166
Definition: dirichletbcenforcement.hh:6
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
auto nonLinearElastic(const MAT &mat)
A helper function to create a non-linear elastic pre finite element.
Definition: nonlinearelastic.hh:330
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
Definition: utils/dirichletvalues.hh:28
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
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:155
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:61
auto getInternalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: nonlinearelastic.hh:172
typename Traits::Basis Basis
Definition: nonlinearelastic.hh:64
typename Traits::Geometry Geometry
Definition: nonlinearelastic.hh:69
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const -> ScalarType
Definition: nonlinearelastic.hh:283
typename Traits::FlatBasis FlatBasis
Definition: nonlinearelastic.hh:65
auto strainFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
The strain function for the given Requirement.
Definition: nonlinearelastic.hh:138
PRE::Material Material
Definition: nonlinearelastic.hh:72
typename Traits::GridView GridView
Definition: nonlinearelastic.hh:70
void bindImpl()
A helper function to bind the local view to the element.
Definition: nonlinearelastic.hh:91
const Geometry & geometry() const
Definition: nonlinearelastic.hh:200
static constexpr int myDim
Definition: nonlinearelastic.hh:77
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Definition: nonlinearelastic.hh:303
NonLinearElastic(const Pre &pre)
Constructor for the NonLinearElastic class.
Definition: nonlinearelastic.hh:84
typename Traits::Element Element
Definition: nonlinearelastic.hh:71
PRE Pre
Definition: nonlinearelastic.hh:73
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:224
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Calculate the matrix associated with the given Requirement.
Definition: nonlinearelastic.hh:257
std::tuple< decltype(makeRT< ResultTypes::PK2Stress >())> SupportedResultTypes
Definition: nonlinearelastic.hh:204
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the material tangent for the given strain.
Definition: nonlinearelastic.hh:154
auto displacementFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Get the displacement function for the given Requirement.
Definition: nonlinearelastic.hh:120
static constexpr auto strainType
Definition: nonlinearelastic.hh:78
typename Traits::LocalView LocalView
Definition: nonlinearelastic.hh:68
int order() const
Definition: nonlinearelastic.hh:202
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: nonlinearelastic.hh:75
auto getStress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: nonlinearelastic.hh:191
size_t numberOfNodes() const
Definition: nonlinearelastic.hh:201
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:42
MAT Material
Definition: nonlinearelastic.hh:43
MAT material
Definition: nonlinearelastic.hh:44
Definition: utils/dirichletvalues.hh:30