version 0.4.1
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>
60class NonLinearElastic : public ResultTypeBase<ResultTypes::PK2Stress>
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 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
156 }
157
166 template <typename ScalarType, int strainDim>
167 auto getInternalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
168 return material<ScalarType>().template storedEnergy<strainType>(strain);
169 }
170
180 template <typename ScalarType, int strainDim, bool voigt = true>
181 auto getStress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
182 return material<ScalarType>().template stresses<strainType, voigt>(strain);
183 }
184
185 const Geometry& geometry() const { return *geo_; }
186 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
187 [[nodiscard]] int order() const { return order_; }
188
189public:
199 template <template <typename, int, int> class RT>
200 requires(supportsResultType<RT>())
202 Dune::PriorityTag<1>) const {
203 using namespace Dune::DerivativeDirections;
204
206 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress>) {
207 const auto uFunction = displacementFunction(req);
208 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
209 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
210
211 return RTWrapper{mat_.template stresses<StrainTags::greenLagrangian>(toVoigt(E))};
212 }
213 }
214
215private:
216 //> CRTP
217 const auto& underlying() const { return static_cast<const FE&>(*this); }
218 auto& underlying() { return static_cast<FE&>(*this); }
219 std::shared_ptr<const Geometry> geo_;
220 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
221 Material mat_;
222 size_t numberOfNodes_{0};
223 int order_{};
224
225 template <typename ScalarType>
226 decltype(auto) material() const {
227 if constexpr (Concepts::AutodiffScalar<ScalarType>)
228 return mat_.template rebind<ScalarType>();
229 else
230 return mat_;
231 }
232
233protected:
241 template <typename ScalarType>
243 const Requirement& par, const MatrixAffordance& affordance, typename Traits::template MatrixType<> K,
244 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
245 using namespace Dune::DerivativeDirections;
246 using namespace Dune;
247 const auto uFunction = displacementFunction(par, dx);
248 const auto eps = strainFunction(par, dx);
249 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
250 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
251 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
252 const auto u = (uFunction.evaluate(gpIndex, on(gridElement))).eval();
253 const auto C = materialTangent(EVoigt);
254
255 const auto stresses = getStress(EVoigt);
256 for (size_t i = 0; i < numberOfNodes_; ++i) {
257 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
258 for (size_t j = 0; j < numberOfNodes_; ++j) {
259 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
260 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
261 K.template block<myDim, myDim>(i * myDim, j * myDim) += (bopI.transpose() * C * bopJ + kgIJ) * intElement;
262 }
263 }
264 }
265 }
266
267 template <typename ScalarType>
269 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx =
270 std::nullopt) const -> ScalarType {
271 using namespace Dune::DerivativeDirections;
272 using namespace Dune;
273
274 const auto eps = strainFunction(par, dx);
275 const auto& lambda = par.parameter();
276 ScalarType energy = 0.0;
277
278 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
279 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
280 const auto internalEnergy = getInternalEnergy(EVoigt);
281 energy += internalEnergy * geo_->integrationElement(gp.position()) * gp.weight();
282 }
283
284 return energy;
285 }
286
287 template <typename ScalarType>
289 const Requirement& par, VectorAffordance affordance, typename Traits::template VectorType<ScalarType> force,
290 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
291 using namespace Dune::DerivativeDirections;
292 using namespace Dune;
293 const auto eps = strainFunction(par, dx);
294
295 // Internal forces
296 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
297 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
298 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
299 const auto stresses = getStress(EVoigt);
300 for (size_t i = 0; i < numberOfNodes_; ++i) {
301 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
302 force.template segment<myDim>(myDim * i) += bopI.transpose() * stresses * intElement;
303 }
304 }
305 }
306};
307
314template <typename MAT>
315auto nonLinearElastic(const MAT& mat) {
317
318 return pre;
319}
320
321} // namespace Ikarus
322
323#else
324 #error NonLinearElastic depends on dune-localfefunctions, which is not included
325#endif
Helper for the autodiff library.
Helper for transform between Dune linear algebra types and Eigen.
Collection of fallback default functions.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
Material property functions and conversion utilities.
Definition of the LinearElastic class for finite element mechanics computations.
Header file for types of loads in Ikarus finite element mechanics.
Definition of several material related enums.
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:166
Definition: assemblermanipulatorbuildingblocks.hh:22
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:315
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: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: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:159
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:272
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:167
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:268
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:185
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:288
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:201
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
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:242
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:187
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:181
size_t numberOfNodes() const
Definition: nonlinearelastic.hh:186
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:32