version 0.4.8
enhancedassumedstrains.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/localfefunctions/derivativetransformators.hh>
14 #include <dune/localfefunctions/linearAlgebraHelper.hh>
15 #include <dune/localfefunctions/meta.hh>
16
25
26namespace Ikarus {
27
28template <typename PreFE, typename FE, typename ES>
29class EnhancedAssumedStrains;
30
35template <typename ES>
37{
39
40 template <typename PreFE, typename FE>
42};
43
55template <typename PreFE, typename FE, typename ESF>
57 : public std::conditional_t<std::same_as<ESF, EAS::LinearStrain>,
58 ResultTypeBase<ResultTypes::linearStress, ResultTypes::linearStressFull>,
59 ResultTypeBase<ResultTypes::PK2Stress, ResultTypes::PK2StressFull,
60 ResultTypes::kirchhoffStress, ResultTypes::cauchyStress>>
61{
62public:
65 using LocalView = typename Traits::LocalView;
66 using Geometry = typename Traits::Geometry;
67 using GridView = typename Traits::GridView;
70
71 static constexpr int myDim = Traits::mydim;
72
73 template <typename ST>
74 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
75
76 template <template <typename, int, int> class RT>
78
83 explicit EnhancedAssumedStrains(const Pre& pre) {
85 "EAS method is only implemented for the total Lagrangian setting.");
86 static_assert(
87 not(FE::strainType == StrainTags::linear) or (std::same_as<EnhancedStrainFunction, EAS::LinearStrain>),
88 "If FE::strainType is linear, then enhancedStrain must also be linear.");
90 }
91
97 bool isDisplacementBased() const { return easVariant_.isDisplacmentBased(); }
98
104 const auto& easVariant() const { return easVariant_; }
105
111 auto numberOfInternalVariables() const { return easVariant_.numberOfInternalVariables(); }
112
127 template <template <typename, int, int> class RT>
128 requires(EnhancedAssumedStrains::template supportsResultType<RT>())
130 Dune::PriorityTag<2>) const {
131 using namespace Dune::DerivativeDirections;
134 const auto ufunc = underlying().displacementFunction(req);
135 auto disp = Dune::viewAsFlatEigenVector(ufunc.coefficientsRef());
136 const auto geo = underlying().localView().element().geometry();
137
138 decltype(auto) mat = [&]() {
139 if constexpr ((isSameResultType<RT, ResultTypes::PK2StressFull> or
140 isSameResultType<RT, ResultTypes::linearStressFull>) and
141 requires { underlying().material().underlying(); })
142 return underlying().material().underlying();
143 else
144 return underlying().material();
145 }();
146
147 RTWrapperType<RT> resultWrapper{};
148 auto calculateAtContribution = [&]<typename EAST>(const EAST& easFunction) {
149 Eigen::VectorXd alpha;
150 alpha.setZero(numberOfInternalVariables());
151 if constexpr (EAST::enhancedStrainSize != 0) {
152 if constexpr (isSameResultType<RT, ResultTypes::linearStressFull> or
153 isSameResultType<RT, ResultTypes::linearStress>) {
154 typename EAST::DType D;
155 calculateDAndLMatrix(easFunction, req, D, L_);
156 alpha = -D.inverse() * L_ * disp;
157 } else
158 alpha = this->alpha_;
159 }
160
161 const auto H = [&]() {
162 if constexpr (std::same_as<EnhancedStrainFunction, EAS::LinearStrain> or
163 std::same_as<EnhancedStrainFunction, EAS::GreenLagrangeStrain>) {
164 return (ufunc.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement))).eval();
165 } else {
166 return EnhancedStrainFunction::computeDisplacementGradient(geo, ufunc, local, easFunction, alpha).eval();
167 }
168 }();
169
170 const auto F = transformStrain<StrainTags::displacementGradient, StrainTags::deformationGradient>(H).eval();
171 const auto enhancedStrain = EnhancedStrainFunction::value(geo, ufunc, local, easFunction, alpha);
172 const auto stress = underlying().calculateStress(mat, enhancedStrain).eval();
173
174 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or isSameResultType<RT, ResultTypes::PK2Stress> or
175 isSameResultType<RT, ResultTypes::linearStressFull> or
176 isSameResultType<RT, ResultTypes::PK2StressFull>)
177 resultWrapper = stress;
178 else if constexpr (isSameResultType<RT, ResultTypes::kirchhoffStress> or
179 isSameResultType<RT, ResultTypes::cauchyStress>) {
180 const auto PK2StressMat = fromVoigt(stress, false).eval();
181 constexpr StressTags toStress =
182 isSameResultType<RT, ResultTypes::kirchhoffStress> ? StressTags::Kirchhoff : StressTags::Cauchy;
183 resultWrapper = toVoigt(transformStress<StressTags::PK2, toStress>(PK2StressMat, F), false);
184 }
185 };
186 easVariant_(calculateAtContribution);
187 return resultWrapper;
188 }
189
199 initializeState();
200 }
201
207 const auto& internalVariable() const { return alpha_; }
208
209protected:
210 void bindImpl() {
211 assert(underlying().localView().bound());
212 easVariant_.bind(underlying().localView().element().geometry());
213 initializeState();
214 }
215
226 const std::remove_reference_t<typename Traits::template VectorType<>>& correction) {
227 using ScalarType = Traits::ctype;
229 auto correctAlpha = [&]<typename EAST>(const EAST& easFunction) {
230 if constexpr (EAST::enhancedStrainSize != 0) {
231 if (correction.size() != par.globalSolution().size())
232 DUNE_THROW(Dune::NotImplemented,
233 "Solution vector and correction vector should be of the same size. Check if DBCOption::Full is "
234 "used. The sizes are " +
235 std::to_string(par.globalSolution().size()) + " and " + std::to_string(correction.size()));
236 const auto& Rtilde = calculateRtilde<ScalarType>(par);
237 const auto localdxBlock = Ikarus::FEHelper::localSolutionBlockVector<Traits, Eigen::VectorXd, double>(
238 correction, underlying().localView());
239 const auto localdx = Dune::viewAsFlatEigenVector(localdxBlock);
240 constexpr int enhancedStrainSize = EAST::enhancedStrainSize;
241 Eigen::Matrix<double, enhancedStrainSize, enhancedStrainSize> D;
242 calculateDAndLMatrix(easFunction, par, D, L_);
243 this->alpha_ -= D.inverse() * (Rtilde + (L_ * localdx));
244 }
245 };
246
247 easVariant_(correctAlpha);
248 }
249
250 inline void easApplicabilityCheck() const {
251 const auto& numberOfNodes = underlying().numberOfNodes();
252 if (not((numberOfNodes == 4 and Traits::mydim == 2) or (numberOfNodes == 9 and Traits::mydim == 2) or
253 (numberOfNodes == 8 and Traits::mydim == 3)) and
254 (not isDisplacementBased()))
255 DUNE_THROW(Dune::NotImplemented, "EAS is only supported for Q1, Q2 and H1 elements");
256 }
257
258 template <typename ScalarType>
259 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
260 typename Traits::template MatrixType<> K,
261 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
262 if (affordance != MatrixAffordance::stiffness)
263 DUNE_THROW(Dune::NotImplemented, "MatrixAffordance not implemented: " + toString(affordance));
265
266 const auto geo = underlying().localView().element().geometry();
267 const auto& uFunction = underlying().displacementFunction(par, dx);
268 const auto& kMFunction = underlying().template materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
269 const auto& kGFunction = underlying().template geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
270 const auto numberOfNodes = underlying().numberOfNodes();
271 const auto& localBasis = underlying().localBasis();
272
273 auto calculateMatrixContribution = [&]<typename EAST>(const EAST& easFunction) {
274 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
275 const auto enhancedStrain = EnhancedStrainFunction::value(geo, uFunction, gp.position(), easFunction, alpha_);
276 const auto stresses = underlying().stress(enhancedStrain);
277 for (size_t i = 0; i < numberOfNodes; ++i) {
278 const auto E_dI = EnhancedStrainFunction::template firstDerivative<0>(geo, uFunction, localBasis, gpIndex,
279 gp.position(), easFunction, alpha_, i);
280 for (size_t j = 0; j < numberOfNodes; ++j) {
281 const auto E_dJ = EnhancedStrainFunction::template firstDerivative<0>(
282 geo, uFunction, localBasis, gpIndex, gp.position(), easFunction, alpha_, j);
283 const auto E_dd = EnhancedStrainFunction::template secondDerivative<0>(
284 geo, uFunction, localBasis, gpIndex, gp.position(), stresses, easFunction, alpha_, i, j);
285 kMFunction(enhancedStrain, E_dI, E_dJ, i, j, gp);
286 kGFunction(E_dd, i, j, gp);
287 }
288 }
289 }
290
291 if constexpr (EAST::enhancedStrainSize != 0) {
292 typename EAST::DType D;
293 calculateDAndLMatrix(easFunction, par, D, L_);
294 K.template triangularView<Eigen::Upper>() -= L_.transpose() * D.inverse() * L_;
295 K.template triangularView<Eigen::StrictlyLower>() = K.transpose();
296 }
297 };
298 easVariant_(calculateMatrixContribution);
299 }
300
301 template <typename ScalarType>
302 inline ScalarType calculateScalarImpl(const Requirement& par, ScalarAffordance affordance,
303 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
306 return underlying().template energyFunction<ScalarType>(par, dx)();
307
308 DUNE_THROW(Dune::NotImplemented,
309 "EAS element do not support any scalar calculations, i.e. they are not derivable from a potential");
310 }
311
312 template <typename ScalarType>
314 typename Traits::template VectorType<ScalarType> force,
315 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
316 if (affordance != VectorAffordance::forces)
317 DUNE_THROW(Dune::NotImplemented, "VectorAffordance not implemented: " + toString(affordance));
319
320 const auto geo = underlying().localView().element().geometry();
321 const auto& uFunction = underlying().displacementFunction(par, dx);
322 const auto numberOfNodes = underlying().numberOfNodes();
323 const auto& localBasis = underlying().localBasis();
324 const auto fIntFunction = underlying().template internalForcesFunction<ScalarType>(par, force, dx);
325
326 auto calculateForceContribution = [&]<typename EAST>(const EAST& easFunction) {
327 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
328 const auto enhancedStrain = EnhancedStrainFunction::value(geo, uFunction, gp.position(), easFunction, alpha_);
329 const auto stresses = underlying().stress(enhancedStrain);
330 for (size_t i = 0; i < numberOfNodes; ++i) {
331 const auto E_dI = EnhancedStrainFunction::template firstDerivative<0>(geo, uFunction, localBasis, gpIndex,
332 gp.position(), easFunction, alpha_, i);
333 fIntFunction(stresses, E_dI, i, gp);
334 }
335 }
336 if constexpr (EAST::enhancedStrainSize != 0) {
337 // Here, L_ and D doesn't see the ScalarType, while Rtilde sees it. This is needed because only then when we use
338 // automatic differentiation w.r.t to the displacements, we get the correct static-condensed stiffness matrix.
339 // If L_ and D sees the ScalarType, then their derivatives w.r.t. the displacements is also performed which is
340 // then not equivalent to the static-condensed stiffness matrix in general.
341 typename EAST::DType D;
342 calculateDAndLMatrix(easFunction, par, D, L_);
343 const auto& Rtilde = calculateRtilde(par, dx);
344 force -= L_.transpose() * D.inverse() * Rtilde;
345 }
346 };
347 easVariant_(calculateForceContribution);
348 }
349
350 template <typename MT, typename BC>
351 void subscribeToImpl(BC& bc) {
352 if constexpr (std::same_as<MT, NonLinearSolverMessages>) {
353 using NLSState = typename BC::State;
354 underlying().subscribe(bc, [&](NonLinearSolverMessages message, const NLSState& state) {
356 this->updateStateImpl(state.domain, state.correction);
357 });
358 }
359 }
360
361private:
363 mutable Eigen::MatrixXd L_;
364 Eigen::VectorXd alpha_;
365
366 //> CRTP
367 const auto& underlying() const { return static_cast<const FE&>(*this); }
368 auto& underlying() { return static_cast<FE&>(*this); }
369
373 void initializeState() {
374 alpha_.resize(numberOfInternalVariables());
375 alpha_.setZero();
376 }
377
378 template <int enhancedStrainSize>
379 void calculateDAndLMatrix(const auto& easFunction, const auto& par,
380 Eigen::Matrix<double, enhancedStrainSize, enhancedStrainSize>& DMat,
381 Eigen::MatrixX<double>& LMat) const {
382 using namespace Dune;
383 using namespace Dune::DerivativeDirections;
384
385 const auto& uFunction = underlying().displacementFunction(par);
386 const auto geo = underlying().localView().element().geometry();
387 const auto numberOfNodes = underlying().numberOfNodes();
388 const auto& localBasis = underlying().localBasis();
389 DMat.setZero();
390 LMat.setZero(enhancedStrainSize, underlying().localView().size());
391 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
392 const auto enhancedStrain = EnhancedStrainFunction::value(geo, uFunction, gp.position(), easFunction, alpha_);
393 const auto C = underlying().materialTangent(enhancedStrain);
394 const auto stresses = underlying().stress(enhancedStrain);
395 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
396 const auto E_a = EnhancedStrainFunction::template firstDerivative<1>(geo, uFunction, localBasis, gpIndex,
397 gp.position(), easFunction, alpha_);
398 const auto E_aa = EnhancedStrainFunction::template secondDerivative<1>(
399 geo, uFunction, localBasis, gpIndex, gp.position(), stresses, easFunction, alpha_);
400 DMat += (E_a.transpose() * C * E_a + E_aa) * intElement;
401 for (size_t i = 0U; i < numberOfNodes; ++i) {
402 const auto E_dI = EnhancedStrainFunction::template firstDerivative<0>(geo, uFunction, localBasis, gpIndex,
403 gp.position(), easFunction, alpha_, i);
404 const auto E_ad = EnhancedStrainFunction::template secondDerivative<2>(
405 geo, uFunction, localBasis, gpIndex, gp.position(), stresses, easFunction, alpha_, i);
406 LMat.template block<enhancedStrainSize, myDim>(0, myDim * i) +=
407 (E_a.transpose() * C * E_dI + E_ad) * intElement;
408 }
409 }
410 }
411
412 template <typename ScalarType>
413 Eigen::VectorX<ScalarType> calculateRtilde(const Requirement& par,
414 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
415 const auto geo = underlying().localView().element().geometry();
416 const auto& uFunction = underlying().displacementFunction(par, dx);
417 const auto& localBasis = underlying().localBasis();
418 Eigen::VectorX<ScalarType> Rtilde;
419 Rtilde.setZero(numberOfInternalVariables());
420
421 auto calculateRtildeContribution = [&]<typename EAST>(const EAST& easFunction) {
422 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
423 const auto enhancedStrain = EnhancedStrainFunction::value(geo, uFunction, gp.position(), easFunction, alpha_);
424 const auto E_a = EnhancedStrainFunction::template firstDerivative<1>(geo, uFunction, localBasis, gpIndex,
425 gp.position(), easFunction, alpha_);
426 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
427 auto stresses = underlying().stress(enhancedStrain);
428 Rtilde += (E_a.transpose() * stresses).eval() * intElement;
429 }
430 };
431
432 easVariant_(calculateRtildeContribution);
433 return Rtilde;
434 }
435};
436
443template <typename ES = EAS::LinearStrain>
444auto eas(int numberOfInternalVariables = 0) {
445 EnhancedAssumedStrainsPre<ES> pre(numberOfInternalVariables);
446
447 return pre;
448}
449
450} // namespace Ikarus
451
452#else
453 #error EnhancedAssumedStrains depends on dune-localfefunctions, which is not included
454#endif
Implementation of transformations for different strain measures.
Implementation of transformations for different stress measures.
Definition of several material related enums.
Definition of the EAS variants.
Definition of the LinearElastic class for finite element mechanics computations.
Enums for observer messages.
StressTags
A strongly typed enum class representing the type of the computed stresses.
Definition: tags.hh:26
auto viewAsFlatEigenVector(Dune::BlockVector< ValueType > &blockedVector)
View Dune::BlockVector as an Eigen::Vector.
Definition: linearalgebrahelper.hh:56
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 eas(int numberOfInternalVariables=0)
A helper function to create an enhanced assumed strain pre finite element.
Definition: enhancedassumedstrains.hh:444
NonLinearSolverMessages
Enum class defining non-linear solver-related messages.
Definition: broadcastermessages.hh:22
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
constexpr std::string toString(DBCOption _e)
Definition: dirichletbcenforcement.hh:8
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
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
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
double ctype
Type used for coordinates.
Definition: fetraits.hh:57
typename Basis::GridView GridView
Type of the grid view.
Definition: fetraits.hh:45
static constexpr int mydim
Dimension of the geometry.
Definition: fetraits.hh:63
Wrapper class for using Enhanced Assumed Strains (EAS) with displacement based elements.
Definition: enhancedassumedstrains.hh:61
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: enhancedassumedstrains.hh:74
FERequirements< FESolutions::displacement, FEParameter::loadfactor > Requirement
Definition: enhancedassumedstrains.hh:64
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: enhancedassumedstrains.hh:259
void easApplicabilityCheck() const
Definition: enhancedassumedstrains.hh:250
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: enhancedassumedstrains.hh:313
auto calculateAtImpl(const Requirement &req, const Dune::FieldVector< double, Traits::mydim > &local, Dune::PriorityTag< 2 >) const
Calculates a requested result at a specific local position using the Enhanced Assumed Strains (EAS) m...
Definition: enhancedassumedstrains.hh:129
void bindImpl()
Definition: enhancedassumedstrains.hh:210
ScalarType calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: enhancedassumedstrains.hh:302
ESF EnhancedStrainFunction
Definition: enhancedassumedstrains.hh:69
bool isDisplacementBased() const
Checks if the element is displacement-based and the EAS is turned off.
Definition: enhancedassumedstrains.hh:97
void setEASType(int numberOfInternalVariables)
Sets the EAS type for 2D elements.
Definition: enhancedassumedstrains.hh:195
EnhancedAssumedStrains(const Pre &pre)
Constructor for Enhanced Assumed Strains elements.
Definition: enhancedassumedstrains.hh:83
typename Traits::Geometry Geometry
Definition: enhancedassumedstrains.hh:66
const auto & internalVariable() const
Gets the internal state variable alpha for the EAS element.
Definition: enhancedassumedstrains.hh:207
static constexpr int myDim
Definition: enhancedassumedstrains.hh:71
void subscribeToImpl(BC &bc)
Definition: enhancedassumedstrains.hh:351
const auto & easVariant() const
Gets the variant representing the type of Enhanced Assumed Strains (EAS).
Definition: enhancedassumedstrains.hh:104
typename Traits::GridView GridView
Definition: enhancedassumedstrains.hh:67
typename Traits::LocalView LocalView
Definition: enhancedassumedstrains.hh:65
void updateStateImpl(const Requirement &par, const std::remove_reference_t< typename Traits::template VectorType<> > &correction)
Updates the internal state variable alpha_ at the end of an iteration before the update of the displa...
Definition: enhancedassumedstrains.hh:225
auto numberOfInternalVariables() const
Gets the number of EAS parameters based on the current EAS type.
Definition: enhancedassumedstrains.hh:111
A PreFE struct for Enhanced Assumed Strains.
Definition: enhancedassumedstrains.hh:37
int numberOfParameters
Definition: enhancedassumedstrains.hh:38
void bind(const GEO &geometry)
Definition: easvariants.hh:95
auto numberOfInternalVariables() const
A helper function to get the number of EAS parameters.
Definition: easvariants.hh:79
bool isDisplacmentBased() const
A helper function to identify if the formulation is purely displacement-based, i.e....
Definition: easvariants.hh:88
void setEASType(int numberOfInternalVariables)
Definition: easvariants.hh:90
Definition: utils/dirichletvalues.hh:38
Concept to check if the underlying strain and stress tag correspond to a total Lagrangian formulation...
Definition: utils/concepts.hh:694
Several concepts.