version 0.4
finiteelements/mechanics/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
43 template <typename Basis_, typename Material_, typename FERequirements_ = FERequirements<>, bool useEigenRef = false>
44 class NonLinearElastic : public PowerBasisFE<typename Basis_::FlatBasis>,
45 public Volume<NonLinearElastic<Basis_, Material_, FERequirements_, useEigenRef>,
46 TraitsFromFE<Basis_, FERequirements_, useEigenRef>>,
47 public Traction<NonLinearElastic<Basis_, Material_, FERequirements_, useEigenRef>,
48 TraitsFromFE<Basis_, FERequirements_, useEigenRef>> {
49 public:
51 using Basis = typename Traits::Basis;
52 using FlatBasis = typename Traits::FlatBasis;
54 using LocalView = typename Traits::LocalView;
55 using Geometry = typename Traits::Geometry;
56 using GridView = typename Traits::GridView;
57 using Element = typename Traits::Element;
59 using BasePowerFE = PowerBasisFE<FlatBasis>; // Handles globalIndices function
60 using Material = Material_;
63 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
64
65 static constexpr int myDim = Traits::mydim;
66 static constexpr auto strainType = StrainTags::greenLagrangian;
67
80 template <typename VolumeLoad = utils::LoadDefault, typename NeumannBoundaryLoad = utils::LoadDefault>
81 NonLinearElastic(const Basis& globalBasis, const typename LocalView::Element& element, const Material& p_mat,
82 VolumeLoad p_volumeLoad = {}, const BoundaryPatch<GridView>* p_neumannBoundary = nullptr,
83 NeumannBoundaryLoad p_neumannBoundaryLoad = {})
84 : BasePowerFE(globalBasis.flat(), element),
85 VolumeType(p_volumeLoad),
86 TractionType(p_neumannBoundary, p_neumannBoundaryLoad),
87 mat{p_mat} {
88 this->localView().bind(element);
89 auto& first_child = this->localView().tree().child(0);
90 const auto& fe = first_child.finiteElement();
91 geo_ = std::make_shared<const Geometry>(this->localView().element().geometry());
92 numberOfNodes_ = fe.size();
93 order_ = 2 * (fe.localBasis().order());
94 localBasis = Dune::CachedLocalBasis(fe.localBasis());
95 if constexpr (requires { this->localView().element().impl().getQuadratureRule(order_); })
96 if (this->localView().element().impl().isTrimmed())
97 localBasis.bind(this->localView().element().impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
98 else
99 localBasis.bind(Dune::QuadratureRules<double, myDim>::rule(this->localView().element().type(), order_),
100 Dune::bindDerivatives(0, 1));
101 else
102 localBasis.bind(Dune::QuadratureRules<double, myDim>::rule(this->localView().element().type(), order_),
103 Dune::bindDerivatives(0, 1));
104 }
105
114 template <typename ScalarType = double>
116 const std::optional<const Eigen::VectorX<ScalarType>>& dx = std::nullopt) const {
117 const auto& d = par.getGlobalSolution(Ikarus::FESolutions::displacement);
118 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, this->localView(), dx);
119 Dune::StandardLocalFunction uFunction(localBasis, disp, geo_);
120 return uFunction;
121 }
122
131 template <typename ScalarType = double>
132 inline auto strainFunction(const FERequirementType& par,
133 const std::optional<const Eigen::VectorX<ScalarType>>& dx = std::nullopt) const {
134 return Dune::greenLagrangeStrains(displacementFunction(par, dx));
135 }
136
146 template <typename ScalarType, int strainDim, bool voigt = true>
147 auto materialTangent(const Eigen::Vector<ScalarType, strainDim>& strain) const {
148 if constexpr (std::is_same_v<ScalarType, double>)
149 return mat.template tangentModuli<strainType, voigt>(strain);
150 else {
151 decltype(auto) matAD = mat.template rebind<ScalarType>();
152 return matAD.template tangentModuli<strainType, voigt>(strain);
153 }
154 }
155
164 template <typename ScalarType, int strainDim>
165 auto getInternalEnergy(const Eigen::Vector<ScalarType, strainDim>& strain) const {
166 if constexpr (std::is_same_v<ScalarType, double>)
167 return mat.template storedEnergy<strainType>(strain);
168 else {
169 decltype(auto) matAD = mat.template rebind<ScalarType>();
170 return matAD.template storedEnergy<strainType>(strain);
171 }
172 }
173
183 template <typename ScalarType, int strainDim, bool voigt = true>
184 auto getStress(const Eigen::Vector<ScalarType, strainDim>& strain) const {
185 if constexpr (std::is_same_v<ScalarType, double>)
186 return mat.template stresses<strainType, voigt>(strain);
187 else {
188 decltype(auto) matAD = mat.template rebind<ScalarType>();
189 return matAD.template stresses<strainType, voigt>(strain);
190 }
191 }
192
193 const Geometry& geometry() const { return *geo_; }
194 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
195 [[nodiscard]] int order() const { return order_; }
196
204 double calculateScalar(const FERequirementType& par) const { return calculateScalarImpl<double>(par); }
205
213 void calculateVector(const FERequirementType& par, typename Traits::template VectorType<> force) const {
214 calculateVectorImpl<double>(par, force);
215 }
216
224 void calculateMatrix(const FERequirementType& par, typename Traits::template MatrixType<> K) const {
225 using namespace Dune::DerivativeDirections;
226 using namespace Dune;
227 const auto eps = strainFunction(par);
228 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
229 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
230 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
231 const auto C = materialTangent(EVoigt);
232 const auto stresses = getStress(EVoigt);
233 for (size_t i = 0; i < numberOfNodes_; ++i) {
234 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
235 for (size_t j = 0; j < numberOfNodes_; ++j) {
236 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
237 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
238 K.template block<myDim, myDim>(i * myDim, j * myDim) += (bopI.transpose() * C * bopJ + kgIJ) * intElement;
239 }
240 }
241 }
242
243 // Update due to displacement-dependent external loads, e.g., follower loads
246 }
247
256 ResultTypeMap<double>& result) const {
257 using namespace Dune::DerivativeDirections;
258 using namespace Dune;
259
260 const auto uFunction = displacementFunction(req.getFERequirements());
261 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
262 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
263 const auto EVoigt = toVoigt(E);
264 auto PK2 = mat.template stresses<StrainTags::greenLagrangian>(EVoigt);
265
266 if (req.isResultRequested(ResultType::PK2Stress))
268 else
269 DUNE_THROW(Dune::NotImplemented, "The requested result type is NOT implemented.");
270 }
271
272 private:
273 std::shared_ptr<const Geometry> geo_;
274 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis;
275 Material mat;
276 size_t numberOfNodes_{0};
277 int order_{};
278
279 protected:
280 template <typename ScalarType>
281 auto calculateScalarImpl(const FERequirementType& par, const std::optional<const Eigen::VectorX<ScalarType>>& dx
282 = std::nullopt) const -> ScalarType {
283 using namespace Dune::DerivativeDirections;
284 using namespace Dune;
285 const auto uFunction = displacementFunction(par, dx);
286 const auto eps = strainFunction(par, dx);
287 const auto& lambda = par.getParameter(Ikarus::FEParameter::loadfactor);
288 ScalarType energy = 0.0;
289
290 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
291 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
292 const auto internalEnergy = getInternalEnergy(EVoigt);
293 energy += internalEnergy * geo_->integrationElement(gp.position()) * gp.weight();
294 }
295
296 // External forces volume forces over the domain
297 energy += VolumeType::calculateScalarImpl(par, dx);
298
299 // line or surface loads, i.e., neumann boundary
300 energy += TractionType::calculateScalarImpl(par, dx);
301 return energy;
302 }
303
304 template <typename ScalarType>
305 void calculateVectorImpl(const FERequirementType& par, typename Traits::template VectorType<ScalarType> force,
306 const std::optional<const Eigen::VectorX<ScalarType>>& dx = std::nullopt) const {
307 using namespace Dune::DerivativeDirections;
308 using namespace Dune;
309 const auto eps = strainFunction(par, dx);
310
311 // Internal forces
312 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
313 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
314 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
315 const auto stresses = getStress(EVoigt);
316 for (size_t i = 0; i < numberOfNodes_; ++i) {
317 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
318 force.template segment<myDim>(myDim * i) += bopI.transpose() * stresses * intElement;
319 }
320 }
321
322 // External forces volume forces over the domain
323 VolumeType::calculateVectorImpl(par, force, dx);
324
325 // External forces, boundary forces, i.e., at the Neumann boundary
326 TractionType::calculateVectorImpl(par, force, dx);
327 }
328 };
329} // namespace Ikarus
330
331#else
332# error NonLinearElastic depends on dune-localfefunctions, which is not included
333#endif
Collection of fallback default functions.
Helper for transform between Dune linear algebra types and Eigen.
Helper for the autodiff library.
Contains the PowerBasisFE class, which works with a power basis in FlatInterLeaved elements.
Material property functions and conversion utilities.
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.
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:167
Definition: simpleassemblers.hh:21
Definition: resultevaluators.hh:17
PowerBasisFE class for working with a power basis in FlatInterLeaved elements.
Definition: powerbasisfe.hh:23
const GridElement & gridElement() const
Get the grid element associated with the local view.
Definition: powerbasisfe.hh:78
const LocalView & localView() const
Get the const reference to the local view.
Definition: powerbasisfe.hh:84
Class representing a map of result types to result arrays.
Definition: ferequirements.hh:342
void insertOrAssignResult(ResultType &&resultType, const ResultArray &resultArray)
Insert or assign a result to the map.
Definition: ferequirements.hh:354
Traction class represents distributed traction load that can be applied.
Definition: traction.hh:22
void calculateVectorImpl(const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > dx=std::nullopt) const
Definition: traction.hh:112
auto calculateScalarImpl(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType
Definition: traction.hh:80
void calculateMatrix(const FERequirementType &req, typename Traits::template MatrixType<> K) const
Calculate the matrix associated with the given FERequirementType.
Definition: traction.hh:74
Volume class represents distributed volume load that can be applied.
Definition: volume.hh:20
auto calculateScalarImpl(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType
Definition: volume.hh:71
void calculateMatrix(const FERequirementType &req, typename Traits::template MatrixType<> K) const
Calculate the matrix associated with the given FERequirementType.
Definition: volume.hh:65
void calculateVectorImpl(const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Definition: volume.hh:88
Interface classf or materials.
Definition: interface.hh:75
NonLinearElastic class represents a non-linear elastic finite element.
Definition: finiteelements/mechanics/nonlinearelastic.hh:48
Traction< NonLinearElastic< Basis_, Material_, FERequirements_, useEigenRef >, Traits > TractionType
Definition: finiteelements/mechanics/nonlinearelastic.hh:62
static constexpr int myDim
Definition: finiteelements/mechanics/nonlinearelastic.hh:65
typename Traits::ResultRequirementsType ResultRequirementsType
Definition: finiteelements/mechanics/nonlinearelastic.hh:58
const Geometry & geometry() const
Definition: finiteelements/mechanics/nonlinearelastic.hh:193
static constexpr auto strainType
Definition: finiteelements/mechanics/nonlinearelastic.hh:66
int order() const
Definition: finiteelements/mechanics/nonlinearelastic.hh:195
PowerBasisFE< FlatBasis > BasePowerFE
Definition: finiteelements/mechanics/nonlinearelastic.hh:59
void calculateMatrix(const FERequirementType &par, typename Traits::template MatrixType<> K) const
Calculate the matrix associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:224
auto getInternalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:165
typename Traits::Element Element
Definition: finiteelements/mechanics/nonlinearelastic.hh:57
Volume< NonLinearElastic< Basis_, Material_, FERequirements_, useEigenRef >, Traits > VolumeType
Definition: finiteelements/mechanics/nonlinearelastic.hh:61
typename Traits::LocalView LocalView
Definition: finiteelements/mechanics/nonlinearelastic.hh:54
double calculateScalar(const FERequirementType &par) const
Calculate the scalar value associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:204
typename Traits::FERequirementType FERequirementType
Definition: finiteelements/mechanics/nonlinearelastic.hh:53
auto displacementFunction(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Get the displacement function for the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:115
auto calculateScalarImpl(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType
Definition: finiteelements/mechanics/nonlinearelastic.hh:281
void calculateVectorImpl(const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Definition: finiteelements/mechanics/nonlinearelastic.hh:305
Material_ Material
Definition: finiteelements/mechanics/nonlinearelastic.hh:60
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: finiteelements/mechanics/nonlinearelastic.hh:63
size_t numberOfNodes() const
Definition: finiteelements/mechanics/nonlinearelastic.hh:194
auto getStress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:184
typename Traits::GridView GridView
Definition: finiteelements/mechanics/nonlinearelastic.hh:56
typename Traits::Geometry Geometry
Definition: finiteelements/mechanics/nonlinearelastic.hh:55
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the material tangent for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:147
auto strainFunction(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
The strain function for the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:132
typename Traits::Basis Basis
Definition: finiteelements/mechanics/nonlinearelastic.hh:51
typename Traits::FlatBasis FlatBasis
Definition: finiteelements/mechanics/nonlinearelastic.hh:52
void calculateAt(const ResultRequirementsType &req, const Dune::FieldVector< double, Traits::mydim > &local, ResultTypeMap< double > &result) const
Calculate specified results at a given local position.
Definition: finiteelements/mechanics/nonlinearelastic.hh:255
NonLinearElastic(const Basis &globalBasis, const typename LocalView::Element &element, const Material &p_mat, VolumeLoad p_volumeLoad={}, const BoundaryPatch< GridView > *p_neumannBoundary=nullptr, NeumannBoundaryLoad p_neumannBoundaryLoad={})
Constructor for the NonLinearElastic class.
Definition: finiteelements/mechanics/nonlinearelastic.hh:81
void calculateVector(const FERequirementType &par, typename Traits::template VectorType<> force) const
Calculate the vector associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:213
Traits for handling finite elements.see https://en.wikipedia.org/wiki/Lam%C3%A9_parameters.
Definition: physicshelper.hh:68
typename LocalView::Element Element
Type of the grid element.
Definition: physicshelper.hh:88
typename FlatBasis::LocalView LocalView
Type of the local view.
Definition: physicshelper.hh:82
typename Basis::FlatBasis FlatBasis
Type of the flat basis.
Definition: physicshelper.hh:79
ResultRequirements< FERequirementType > ResultRequirementsType
Type of the result requirements.
Definition: physicshelper.hh:76
Basis_ Basis
Type of the basis of the finite element.
Definition: physicshelper.hh:70
static constexpr int mydim
Dimension of the geometry.
Definition: physicshelper.hh:97
typename FlatBasis::GridView GridView
Type of the grid view.
Definition: physicshelper.hh:85
FERequirements_ FERequirementType
Type of the requirements for the finite element.
Definition: physicshelper.hh:73
typename Element::Geometry Geometry
Type of the element geometry.
Definition: physicshelper.hh:91
Definition: resultevaluators.hh:20