version 0.4.8
assumedstress.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
27
28namespace Ikarus {
29
30template <typename PreFE, typename FE, typename ASType>
31class AssumedStress;
32
37template <typename ASType>
39{
41
42 template <typename PreFE, typename FE>
44};
45
57template <typename PreFE, typename FE, typename ASF>
59 : public std::conditional_t<std::same_as<ASF, PS::LinearStress>,
60 ResultTypeBase<ResultTypes::linearStress, ResultTypes::linearStressFull>,
61 ResultTypeBase<ResultTypes::PK2Stress, ResultTypes::PK2StressFull,
62 ResultTypes::kirchhoffStress, ResultTypes::cauchyStress>>
63
64{
65public:
68 using LocalView = typename Traits::LocalView;
69 using Geometry = typename Traits::Geometry;
70 using GridView = typename Traits::GridView;
73
74 template <typename ST>
75 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
76
77 template <template <typename, int, int> class RT>
79
80 static constexpr int myDim = Traits::mydim;
81 static constexpr int strainDim = myDim * (myDim + 1) / 2;
82 using StrainVector = Eigen::Vector<double, strainDim>;
83 using MaterialMatrix = Eigen::Matrix<double, strainDim, strainDim>;
84
89 explicit AssumedStress(const Pre& pre) {
91 "PS method is only implemented for the total Lagrangian setting.");
92 static_assert(not(FE::strainType == StrainTags::linear) or (std::same_as<ASF, PS::LinearStress>),
93 "If FE::strainType is linear, then the assumed stress must also be linear.");
95 }
96
102 const auto& asVariant() const { return asVariant_; }
103
109 auto numberOfInternalVariables() const { return asVariant_.numberOfInternalVariables(); }
110
124 template <template <typename, int, int> class RT>
125 requires(AssumedStress::template supportsResultType<RT>())
127 Dune::PriorityTag<2>) const {
128 using namespace Dune::DerivativeDirections;
131 const auto geo = underlying().localView().element().geometry();
132 const auto ufunc = underlying().displacementFunction(req);
133 auto disp = Dune::viewAsFlatEigenVector(ufunc.coefficientsRef());
134 const auto H = ufunc.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
135 const auto F = transformStrain<StrainTags::displacementGradient, StrainTags::deformationGradient>(H).eval();
136
137 RTWrapperType<RT> resultWrapper{};
138 auto calculateAtContribution = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
139 Eigen::VectorXd beta;
140 beta.setZero(numberOfInternalVariables());
141 if constexpr (isSameResultType<RT, ResultTypes::linearStressFull> or
142 isSameResultType<RT, ResultTypes::linearStress>) {
143 const auto& Rtilde = calculateRtilde<double>(req);
144 typename AssumedStressT::HType H;
145 calculateHAndGMatrix(asFunction, req, H, G_);
146 beta = H.inverse() * (G_ * disp);
147 } else
148 beta = this->beta_;
149
150 const auto SVoigt = AssumedStressFunction::value(geo, local, asFunction, beta);
151
152 if constexpr ((isSameResultType<RT, ResultTypes::PK2StressFull> or
153 isSameResultType<RT, ResultTypes::linearStressFull>) and
154 requires { underlying().material().underlying(); }) {
155 resultWrapper = enlargeStressAnsatz(SVoigt);
156 } else if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or
157 isSameResultType<RT, ResultTypes::linearStress>)
158 resultWrapper = SVoigt;
159 else if constexpr (isSameResultType<RT, ResultTypes::kirchhoffStress> or
160 isSameResultType<RT, ResultTypes::cauchyStress>) {
161 const auto PK2StressMat = fromVoigt(SVoigt, false).eval();
162 constexpr StressTags toStress =
163 isSameResultType<RT, ResultTypes::kirchhoffStress> ? StressTags::Kirchhoff : StressTags::Cauchy;
164 resultWrapper = toVoigt(transformStress<StressTags::PK2, toStress>(PK2StressMat, F), false);
165 }
166 };
167 asVariant_(calculateAtContribution);
168 return resultWrapper;
169 }
170
179 initializeState();
180 }
181
187 const auto& internalVariable() const { return beta_; }
188
189protected:
190 void bindImpl() {
191 assert(underlying().localView().bound());
192 asVariant_.bind(underlying().localView().element().geometry());
193 initializeState();
194 }
195
206 const std::remove_reference_t<typename Traits::template VectorType<>>& correction) {
207 using ScalarType = Traits::ctype;
209 auto correctbeta = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
210 if (correction.size() != par.globalSolution().size())
211 DUNE_THROW(Dune::NotImplemented,
212 "Solution vector and correction vector should be of the same size. Check if DBCOption::Full is "
213 "used. The sizes are " +
214 std::to_string(par.globalSolution().size()) + " and " + std::to_string(correction.size()));
215 const auto localdxBlock = Ikarus::FEHelper::localSolutionBlockVector<Traits, Eigen::VectorXd, double>(
216 correction, underlying().localView());
217 const auto localdx = Dune::viewAsFlatEigenVector(localdxBlock);
218
219 typename AssumedStressT::HType H;
220 calculateHAndGMatrix(asFunction, par, H, G_);
221 const auto& Rtilde = calculateRtilde<ScalarType>(par);
222 this->beta_ += H.inverse() * (Rtilde + (G_ * localdx));
223 };
224
225 asVariant_(correctbeta);
226 calculateMaterialInversion();
227 }
228
229 inline void asApplicabilityCheck() const {
230 const auto& numberOfNodes = underlying().numberOfNodes();
231 if (not((numberOfNodes == 4 and Traits::mydim == 2) or (numberOfNodes == 8 and Traits::mydim == 3)))
232 DUNE_THROW(Dune::NotImplemented, "AssumedStress is only supported for Q1 and H1 elements" +
233 std::to_string(numberOfNodes) + std::to_string(Traits::mydim));
234 }
235
236 template <typename ScalarType>
237 void calculateMatrixImpl(const Requirement& par, const MatrixAffordance& affordance,
238 typename Traits::template MatrixType<> K,
239 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
240 using namespace Dune;
241 using namespace Dune::DerivativeDirections;
242
243 if (affordance != MatrixAffordance::stiffness)
244 DUNE_THROW(Dune::NotImplemented, "MatrixAffordance not implemented: " + toString(affordance));
246
247 const auto geo = underlying().localView().element().geometry();
248 const auto& strainFunction = underlying().strainFunction(par, dx);
249 const auto& kGFunction = underlying().template geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
250 const auto numberOfNodes = underlying().numberOfNodes();
251 const auto& localBasis = underlying().localBasis();
252
253 auto calculateMatrixContribution = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
254 for (const auto& [gpIndex, gp] : strainFunction.viewOverIntegrationPoints()) {
255 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
256
257 for (size_t i = 0; i < numberOfNodes; ++i)
258 for (size_t j = 0; j < numberOfNodes; ++j) {
259 const auto E_dd =
260 strainFunction.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(SVoigt), on(gridElement));
261 kGFunction(E_dd, i, j, gp);
262 }
263 }
264
265 typename AssumedStressT::HType H;
266 calculateHAndGMatrix(asFunction, par, H, G_, dx);
267
268 K.template triangularView<Eigen::Upper>() += G_.transpose() * H.inverse() * G_;
269 K.template triangularView<Eigen::StrictlyLower>() = K.transpose();
270 };
271 asVariant_(calculateMatrixContribution);
272 }
273
274 template <typename ScalarType>
275 inline ScalarType calculateScalarImpl(const Requirement& par, ScalarAffordance affordance,
276 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
277 using namespace Dune;
278 using namespace Dune::DerivativeDirections;
280 ScalarType energy = 0.0;
281
282 // This will not work in the context of autodiff because `AutoDiffFE` class doesn't include static condensation of
283 // the internal variable
284 if constexpr (not Concepts::AutodiffScalar<ScalarType>) {
285 const auto geo = underlying().localView().element().geometry();
286 const auto strainFunction = underlying().strainFunction(par, dx);
287
288 auto calculateScalarContribution = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
289 for (const auto& [gpIndex, gp] : strainFunction.viewOverIntegrationPoints()) {
290 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
291 const auto& Es = strainStates_[gpIndex];
292 const auto EVoigt = strainFunction.evaluate(gp.position(), on(gridElement));
293
294 energy += ((SVoigt.transpose() * EVoigt - ScalarType(0.5) * (SVoigt.transpose() * Es)).eval() *
295 geo.integrationElement(gp.position()) * gp.weight())
296 .eval()[0];
297 }
298 };
299
300 asVariant_(calculateScalarContribution);
301 } else {
302 DUNE_THROW(Dune::NotImplemented,
303 "AssumedStress element does not support scalar calculations for autodiff scalars");
304 }
305 return energy;
306 }
307
308 template <typename ScalarType>
310 typename Traits::template VectorType<ScalarType> force,
311 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
312 using namespace Dune::DerivativeDirections;
313 using namespace Dune;
314
315 if (affordance != VectorAffordance::forces)
316 DUNE_THROW(Dune::NotImplemented, "VectorAffordance not implemented: " + toString(affordance));
318
319 const auto geo = underlying().localView().element().geometry();
320 const auto& uFunction = underlying().displacementFunction(par, dx);
321 const auto& strainFunction = underlying().strainFunction(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 AssumedStressT>(const AssumedStressT& asFunction) {
327 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
328 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
329
330 for (size_t i = 0; i < numberOfNodes; ++i) {
331 const auto E_dI = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
332 fIntFunction(SVoigt, E_dI, i, gp);
333 }
334 }
335 typename AssumedStressT::HType H;
336 calculateHAndGMatrix(asFunction, par, H, G_);
337 const auto Rtilde = calculateRtilde(par, dx);
338 force += G_.transpose() * H.inverse() * Rtilde;
339 };
340 asVariant_(calculateForceContribution);
341 }
342
343 template <typename MT, typename BC>
344 void subscribeToImpl(BC& bc) {
345 if constexpr (std::same_as<MT, NonLinearSolverMessages>) {
346 using NLSState = typename BC::State;
347 underlying().subscribe(bc, [&](NonLinearSolverMessages message, const NLSState& state) {
349 this->updateStateImpl(state.domain, state.correction);
350 });
351 }
352 }
353
354private:
356 mutable Eigen::MatrixXd G_;
357 Eigen::VectorXd beta_;
358 std::vector<StrainVector> strainStates_;
359 std::vector<MaterialMatrix> invertedMaterialStates_;
360
361 //> CRTP
362 const auto& underlying() const { return static_cast<const FE&>(*this); }
363 auto& underlying() { return static_cast<FE&>(*this); }
364
368 void initializeState() {
369 beta_.resize(numberOfInternalVariables());
370 beta_.setZero();
371 strainStates_.resize(underlying().localBasis().integrationPointSize(), StrainVector::Zero());
372 invertedMaterialStates_.resize(underlying().localBasis().integrationPointSize(), MaterialMatrix::Zero());
373 calculateMaterialInversion();
374 }
375
376 template <typename ScalarType, int assumedStressSize>
377 void calculateHAndGMatrix(const auto& asFunction, const auto& par,
378 Eigen::Matrix<double, assumedStressSize, assumedStressSize>& H,
379 Eigen::MatrixX<ScalarType>& G, const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
380 using namespace Dune::DerivativeDirections;
381 using namespace Dune;
382
383 const auto& uFunction = underlying().displacementFunction(par);
384 const auto& strainFunction = underlying().strainFunction(par, dx);
385 const auto geo = underlying().localView().element().geometry();
386 const auto numberOfNodes = underlying().numberOfNodes();
387 const auto& localBasis = underlying().localBasis();
388 H.setZero();
389 G.setZero(assumedStressSize, underlying().localView().size());
390 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
391 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
392
393 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
394 const auto& D = invertedMaterialStates_[gpIndex];
395
396 const auto S_b = AssumedStressFunction::template firstDerivative(geo, uFunction, localBasis, gpIndex,
397 gp.position(), asFunction, beta_);
398 const auto S_bb = AssumedStressFunction::template secondDerivative(geo, uFunction, localBasis, gpIndex,
399 gp.position(), SVoigt, asFunction, beta_);
400
401 H += (S_b.transpose() * D * S_b + S_bb) * intElement;
402
403 for (size_t i = 0U; i < numberOfNodes; ++i) {
404 const auto E_dI = strainFunction.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
405
406 G.template block<assumedStressSize, Traits::worlddim>(0, myDim * i) += S_b.transpose() * E_dI * intElement;
407 }
408 }
409 }
410
411 template <typename ScalarType>
412 Eigen::VectorX<ScalarType> calculateRtilde(const Requirement& par,
413 const VectorXOptRef<ScalarType>& dx = std::nullopt) const {
414 using namespace Dune;
415 using namespace Dune::DerivativeDirections;
416
417 const auto geo = underlying().localView().element().geometry();
418 const auto& uFunction = underlying().displacementFunction(par, dx);
419 const auto strainFunction = underlying().strainFunction(par, dx);
420 const auto& localBasis = underlying().localBasis();
421
422 Eigen::VectorX<ScalarType> Rtilde;
423 Rtilde.setZero(numberOfInternalVariables());
424
425 auto calculateRtildeContribution = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
426 for (const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
427 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
428 const auto S_b = AssumedStressFunction::template firstDerivative(geo, uFunction, localBasis, gpIndex,
429 gp.position(), asFunction, beta_);
430
431 const auto& Es = strainStates_[gpIndex];
432 const auto EVoigt = strainFunction.evaluate(gp.position(), on(gridElement));
433
434 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
435 Rtilde += (S_b.transpose() * (EVoigt - Es)).eval() * intElement;
436 }
437 };
438
439 asVariant_(calculateRtildeContribution);
440 return Rtilde;
441 }
442
443 void calculateMaterialInversion() {
444 auto saveMaterialState = [&]<typename AssumedStressT>(const AssumedStressT& asFunction) {
445 const auto geo = underlying().localView().element().geometry();
446
447 for (const auto& [gpIndex, gp] : underlying().localBasis().viewOverIntegrationPoints()) {
448 const auto SVoigt = AssumedStressFunction::value(geo, gp.position(), asFunction, beta_);
449 const auto& EsOld = strainStates_[gpIndex];
450 std::tie(invertedMaterialStates_[gpIndex], strainStates_[gpIndex]) =
451 underlying().material().template materialInversion<FE::strainType, true>(SVoigt, EsOld);
452 }
453 };
454 asVariant_(saveMaterialState);
455 }
456
457 auto enlargeStressAnsatz(const StrainVector& SVoigt, size_t qpIndexForIteration = 0) const
458 requires(decltype(underlying().material())::isReduced)
459 {
460 auto [_, Es] = underlying().material().template materialInversion<FE::strainType, true>(
461 SVoigt, strainStates_[qpIndexForIteration]);
462 auto Esfull = enlargeIfReduced<decltype(underlying().material())>(Es).eval();
463 return underlying().material().underlying().template stresses<FE::strainType, true>(Esfull);
464 }
465};
466
472template <typename ASType = PS::LinearStress>
473auto assumedStress(int numberOfInternalVariables) {
474 return AssumedStressPre<ASType>(numberOfInternalVariables);
475}
476
477} // namespace Ikarus
478
479#else
480 #error AssumedStress depends on dune-localfefunctions, which is not included
481#endif
Helper for the autodiff library.
Header file for various assumed stress functions.
Definition of the AssumedStress variants.
Implementation of transformations for different strain measures.
Implementation of transformations for different stress measures.
Definition of several material related enums.
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
decltype(auto) enlargeIfReduced(const Eigen::MatrixBase< Derived > &E)
Enlarges a matrix if it reduced in the context of material laws, i.e., VanishingStress If the materia...
Definition: linearalgebrahelper.hh:617
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
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
auto assumedStress(int numberOfInternalVariables)
A helper function to create an assumed stress pre finite element.
Definition: assumedstress.hh:473
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 Assumed Stress with displacement based elements.
Definition: assumedstress.hh:64
const auto & internalVariable() const
Gets the internal state variable beta for the AssumedStress element.
Definition: assumedstress.hh:187
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: assumedstress.hh:237
void subscribeToImpl(BC &bc)
Definition: assumedstress.hh:344
typename Traits::LocalView LocalView
Definition: assumedstress.hh:68
AssumedStress(const Pre &pre)
Constructor for Assunmed Stress elements.
Definition: assumedstress.hh:89
static constexpr int myDim
Definition: assumedstress.hh:80
void setAssumedStressType(int numberOfInternalVariables)
Sets the AssumedStress type for 2D elements.
Definition: assumedstress.hh:176
typename Traits::Geometry Geometry
Definition: assumedstress.hh:69
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 Assumed Stress Method.
Definition: assumedstress.hh:126
void updateStateImpl(const Requirement &par, const std::remove_reference_t< typename Traits::template VectorType<> > &correction)
Updates the internal state variable beta_ at the end of an iteration before the update of the displac...
Definition: assumedstress.hh:205
typename Traits::GridView GridView
Definition: assumedstress.hh:70
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: assumedstress.hh:309
Eigen::Matrix< double, strainDim, strainDim > MaterialMatrix
Definition: assumedstress.hh:83
FERequirements< FESolutions::displacement, FEParameter::loadfactor > Requirement
Definition: assumedstress.hh:67
void bindImpl()
Definition: assumedstress.hh:190
ASF AssumedStressFunction
Definition: assumedstress.hh:72
ScalarType calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: assumedstress.hh:275
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: assumedstress.hh:75
static constexpr int strainDim
Definition: assumedstress.hh:81
auto numberOfInternalVariables() const
Gets the number of AssumedStress parameters based on the current AssumedStress type.
Definition: assumedstress.hh:109
const auto & asVariant() const
Gets the variant representing the type of Assumed Stress.
Definition: assumedstress.hh:102
Eigen::Vector< double, strainDim > StrainVector
Definition: assumedstress.hh:82
void asApplicabilityCheck() const
Definition: assumedstress.hh:229
A PreFE struct for Assumed Stress.
Definition: assumedstress.hh:39
int numberOfParameters
Definition: assumedstress.hh:40
auto numberOfInternalVariables() const
A helper function to get the number of AssumedStress parameters.
Definition: asvariants.hh:63
void bind(const GEO &geometry)
Definition: asvariants.hh:73
void setAssumedStressType(int numberOfInternalVariables)
Definition: asvariants.hh:68
Definition: utils/dirichletvalues.hh:38
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:625
Concept to check if the underlying strain and stress tag correspond to a total Lagrangian formulation...
Definition: utils/concepts.hh:694
Several concepts.