version 0.4.2
linearelastic.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
14 #include <optional>
15 #include <type_traits>
16
17 #include <dune/geometry/quadraturerules.hh>
18 #include <dune/localfefunctions/expressions/linearStrainsExpr.hh>
19 #include <dune/localfefunctions/impl/standardLocalFunction.hh>
20
25
26namespace Ikarus {
27
28template <typename PreFE, typename FE>
29class LinearElastic;
30
35{
37
38 template <typename PreFE, typename FE>
40};
41
50template <typename PreFE, typename FE>
52{
53public:
56 using FlatBasis = typename Traits::FlatBasis;
59 using LocalView = typename Traits::LocalView;
60 using Geometry = typename Traits::Geometry;
61 using GridView = typename Traits::GridView;
62 using Element = typename Traits::Element;
64
65 static constexpr int myDim = Traits::mydim;
66 using LocalBasisType = decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
67
72 explicit LinearElastic(const Pre& pre)
73 : mat_{pre.material} {}
74
75protected:
79 void bindImpl() {
80 const auto& localView = underlying().localView();
81 const auto& element = localView.element();
82 auto& firstChild = localView.tree().child(0);
83 const auto& fe = firstChild.finiteElement();
84 geo_ = std::make_shared<const Geometry>(element.geometry());
85 numberOfNodes_ = fe.size();
86 order_ = 2 * (fe.localBasis().order());
87 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
88 if constexpr (requires { element.impl().getQuadratureRule(order_); })
89 if (element.impl().isTrimmed())
90 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
91 else
92 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
93 Dune::bindDerivatives(0, 1));
94 else
95 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
96 }
97
98public:
107 template <typename ScalarType = double>
109 const Requirement& par,
110 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
111 const auto& d = par.globalSolution();
112 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
113 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
114 return uFunction;
115 }
124 template <class ScalarType = double>
126 const Requirement& par,
127 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
128 return Dune::linearStrains(displacementFunction(par, dx));
129 }
130
136 auto materialTangent() const {
137 if constexpr (myDim == 2)
139 else if constexpr (myDim == 3)
140 return linearElasticMaterialTangent3D(mat_.emodul, mat_.nu);
141 }
142
149 auto materialTangentFunction([[maybe_unused]] const Requirement& par) const {
150 return [&]([[maybe_unused]] auto gp) { return materialTangent(); };
151 }
152
153 const Geometry& geometry() const { return *geo_; }
154 [[nodiscard]] size_t numberOfNodes() const { return numberOfNodes_; }
155 [[nodiscard]] int order() const { return order_; }
156
157 using SupportedResultTypes = std::tuple<decltype(makeRT<ResultTypes::linearStress>())>;
158
159 template <template <typename, int, int> class RT>
160 static consteval bool canProvideResultType() {
161 return isSameResultType<RT, ResultTypes::linearStress>;
162 }
163
164public:
175 template <template <typename, int, int> class RT>
176 requires(canProvideResultType<RT>())
178 Dune::PriorityTag<1>) const {
180 if constexpr (isSameResultType<RT, ResultTypes::linearStress>) {
181 const auto eps = strainFunction(req);
182 const auto C = materialTangent();
183 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
184
185 return RTWrapper{(C * epsVoigt).eval()};
186 }
187 }
188
189private:
190 //> CRTP
191 const auto& underlying() const { return static_cast<const FE&>(*this); }
192 auto& underlying() { return static_cast<FE&>(*this); }
193
194 std::shared_ptr<const Geometry> geo_;
195 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
196 YoungsModulusAndPoissonsRatio mat_;
197 size_t numberOfNodes_{0};
198 int order_{};
199
200protected:
201 template <typename ScalarType>
203 const Requirement& par, const MatrixAffordance& affordance, typename Traits::template MatrixType<> K,
204 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
205 const auto eps = strainFunction(par, dx);
206 using namespace Dune::DerivativeDirections;
207 using namespace Dune;
208
209 const auto C = materialTangent();
210 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
211 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
212 for (size_t i = 0; i < numberOfNodes_; ++i) {
213 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
214 for (size_t j = 0; j < numberOfNodes_; ++j) {
215 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
216 K.template block<myDim, myDim>(i * myDim, j * myDim) += bopI.transpose() * C * bopJ * intElement;
217 }
218 }
219 }
220 }
221
222 template <typename ScalarType>
224 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx =
225 std::nullopt) const -> ScalarType {
226 const auto uFunction = displacementFunction(par, dx);
227 const auto eps = strainFunction(par, dx);
228 const auto& lambda = par.parameter();
229 using namespace Dune::DerivativeDirections;
230 using namespace Dune;
231
232 const auto C = materialTangent();
233
234 ScalarType energy = 0.0;
235 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
236 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
237 energy += (0.5 * EVoigt.dot(C * EVoigt)) * geo_->integrationElement(gp.position()) * gp.weight();
238 }
239 return energy;
240 }
241
242 template <typename ScalarType>
244 const Requirement& par, VectorAffordance affordance, typename Traits::template VectorType<ScalarType> force,
245 const std::optional<std::reference_wrapper<const Eigen::VectorX<ScalarType>>>& dx = std::nullopt) const {
246 const auto eps = strainFunction(par, dx);
247 using namespace Dune::DerivativeDirections;
248 using namespace Dune;
249
250 const auto C = materialTangent();
251
252 // Internal forces
253 for (const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
254 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
255 auto stresses = (C * eps.evaluate(gpIndex, on(gridElement))).eval();
256 for (size_t i = 0; i < numberOfNodes_; ++i) {
257 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
258 force.template segment<myDim>(myDim * i) += bopI.transpose() * stresses * intElement;
259 }
260 }
261 }
262};
263
270 LinearElasticPre pre(mat);
271
272 return pre;
273}
274} // namespace Ikarus
275
276#else
277 #error LinearElastic depends on dune-localfefunctions, which is not included
278#endif
Definition of the LinearElastic class for finite element mechanics computations.
Header file for material models in Ikarus finite element mechanics.
Material property functions and conversion utilities.
Definition: dirichletbcenforcement.hh:6
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
Eigen::Matrix3d planeStressLinearElasticMaterialTangent(double E, double nu)
Computes the plane stress linear elastic material tangent matrix.
Definition: physicshelper.hh:23
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
Eigen::Matrix< double, 6, 6 > linearElasticMaterialTangent3D(double E, double nu)
Computes the 3D linear elastic material tangent matrix.
Definition: physicshelper.hh:40
auto linearElastic(const YoungsModulusAndPoissonsRatio &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:269
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
BH BasisHandler
Type of the basis of the finite element.
Definition: fetraits.hh:27
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
LinearElastic class represents a linear elastic finite element.
Definition: linearelastic.hh:52
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:56
typename Traits::Geometry Geometry
Definition: linearelastic.hh:60
typename Traits::Element Element
Definition: linearelastic.hh:62
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const -> ScalarType
Definition: linearelastic.hh:223
auto materialTangent() const
Gets the material tangent matrix for the linear elastic material.
Definition: linearelastic.hh:136
auto materialTangentFunction(const Requirement &par) const
Gets the material tangent function for the given Requirement.
Definition: linearelastic.hh:149
int order() const
Definition: linearelastic.hh:155
static consteval bool canProvideResultType()
Definition: linearelastic.hh:160
typename Traits::LocalView LocalView
Definition: linearelastic.hh:59
typename Traits::GridView GridView
Definition: linearelastic.hh:61
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:72
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: linearelastic.hh:243
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:55
auto strainFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Gets the strain function for the given Requirement and optional displacement vector.
Definition: linearelastic.hh:125
const Geometry & geometry() const
Definition: linearelastic.hh:153
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:66
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
Definition: linearelastic.hh:202
std::tuple< decltype(makeRT< ResultTypes::linearStress >())> SupportedResultTypes
Definition: linearelastic.hh:157
void bindImpl()
A helper function to bind the local view to the element.
Definition: linearelastic.hh:79
static constexpr int myDim
Definition: linearelastic.hh:65
size_t numberOfNodes() const
Definition: linearelastic.hh:154
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: linearelastic.hh:177
auto displacementFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Gets the displacement function for the given Requirement and optional displacement vector.
Definition: linearelastic.hh:108
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:35
YoungsModulusAndPoissonsRatio material
Definition: linearelastic.hh:36
Structure representing Young's modulus and shear modulus.
Definition: physicshelper.hh:53
double emodul
Definition: physicshelper.hh:54
double nu
Definition: physicshelper.hh:55
Definition: utils/dirichletvalues.hh:30