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>
48 template <Concepts::Material DEV, Concepts::Material VOL>
49 struct DPConstitutiveDriver
53 auto Uhat(
const auto& Estrain)
const {
return vol_.template storedEnergy<strainType>(Estrain); }
55 auto storedEnergy(
const auto& Estrain,
const auto& p) {
56 auto Wdev = dev_.template storedEnergy<strainType>(Estrain);
57 auto uhat = Uhat(Estrain);
58 return Wdev + p * uhat - (1.0 / (2.0 * kappa_)) * Dune::power(p, 2);
61 auto firstDerivatives(
const auto& Estrain,
const auto& p)
const {
62 const auto Sdev = dev_.template stresses<strainType, true>(Estrain);
63 const auto dUdE = vol_.template stresses<strainType, true>(Estrain);
64 const auto S = (Sdev + p * dUdE).eval();
65 return std::make_tuple(dUdE, S);
68 auto secondDerivative(
const auto& Estrain,
const auto& p)
const {
69 const auto Cdev = dev_.template tangentModuli<strainType, true>(Estrain);
70 const auto d2UdE2 = vol_.template tangentModuli<strainType, true>(Estrain);
71 const auto C = (Cdev + p * d2UdE2).eval();
75 double kappa()
const {
return kappa_; }
77 DEV deviatoricMaterial()
const {
return dev_; }
78 VOL volumetricMaterial()
const {
return vol_; }
80 DPConstitutiveDriver(
const DEV& dev,
const VOL& vol,
double kappa)
90 template <
typename STO>
92 auto reboundDEV = dev_.template rebind<STO>();
93 auto reboundVOL = vol_.template rebind<STO>();
94 return DPConstitutiveDriver<decltype(reboundDEV), decltype(reboundVOL)>(reboundDEV, reboundVOL, kappa_);
104template <
typename PreFE,
typename FE,
typename PRE>
105class DisplacementPressure;
111template <Concepts::Material MAT>
117 template <
typename PreFE,
typename FE>
130template <
typename PreFE,
typename FE,
typename PRE>
132 ResultTypes::kirchhoffStress, ResultTypes::cauchyStress>
146 decltype(std::declval<LocalView>().tree().child(Dune::Indices::_0).child(0).finiteElement().localBasis());
148 decltype(std::declval<LocalView>().tree().child(Dune::Indices::_1).finiteElement().localBasis());
150 template <
typename ST>
151 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
158 template <
template <
typename,
int,
int>
class RT>
161 template <
typename ST =
double>
164 template <
typename ST =
double>
165 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
167 template <
typename ST =
double>
168 using KgType = Eigen::Matrix<ST, myDim, myDim>;
172 template <
typename MAT>
175 template <
typename MAT>
186 : constitutiveDriver_([](const auto& mat) {
190 static_assert(Pre::Material::isHyperelastic,
"DisplacementPressure is only implemented for the hyperelastic case.");
198 const auto& localView = underlying().localView();
199 const auto& element = localView.element();
200 auto& firstChild = localView.tree().child(Dune::Indices::_0);
201 const auto& firstFE = firstChild.child(0).finiteElement();
203 auto& secondChild = localView.tree().child(Dune::Indices::_1);
204 const auto& secondFE = secondChild.finiteElement();
206 geo_ = std::make_shared<const Geometry>(element.geometry());
207 numberOfNodes_ = firstFE.size();
208 order_ = 2 * firstFE.localBasis().order();
209 localBasisU_ = Dune::CachedLocalBasis(firstFE.localBasis());
210 localBasisP_ = Dune::CachedLocalBasis(secondFE.localBasis());
212 if (underlying().quadratureRule()) {
213 localBasisU_.bind(underlying().quadratureRule().value(), Dune::bindDerivatives(0, 1));
214 localBasisP_.bind(underlying().quadratureRule().value(), Dune::bindDerivatives(0, 1));
216 localBasisU_.bind(defaultQuadratureRule<myDim>(element, order_), Dune::bindDerivatives(0, 1));
217 localBasisP_.bind(defaultQuadratureRule<myDim>(element, order_), Dune::bindDerivatives(0, 1));
230 template <
typename ScalarType =
double>
234 Ikarus::FEHelper::localSolutionBlockVectorComposite<Traits>(d, underlying().localView(), Dune::Indices::_0, dx);
235 Dune::StandardLocalFunction uFunction(localBasisU_, disp, geo_);
247 template <
typename ScalarType =
double>
251 Ikarus::FEHelper::localSolutionBlockVectorScalar<Traits>(d, underlying().localView(), Dune::Indices::_1, dx);
253 Dune::StandardLocalFunction pFunction(localBasisP_, pressure, geo_);
265 template <
typename ScalarType =
double>
272 [[nodiscard]]
int order()
const {
return order_; }
273 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisTypeU>>&
localBasis()
const {
return localBasisU_; }
274 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisTypeP>>&
localBasisP()
const {
return localBasisP_; }
276 template <
typename ScalarType =
double>
279 return constitutiveDriver_.template rebind<ScalarType>();
281 return constitutiveDriver_;
299 template <
typename D,
typename V>
300 auto calculateStress(
const D& dev,
const V& vol,
const Eigen::Vector<double, strainDim>& strainInVoigt,
302 return (dev.template stresses<strainType>(enlargeIfReduced<typename PRE::Material>(strainInVoigt)) +
303 p * vol.template stresses<strainType>(enlargeIfReduced<typename PRE::Material>(strainInVoigt)))
316 template <
template <
typename,
int,
int>
class RT>
317 requires(supportsResultType<RT>())
319 Dune::PriorityTag<1>)
const {
320 using namespace Dune::DerivativeDirections;
324 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
325 const auto F = transformStrain<StrainTags::displacementGradient, StrainTags::deformationGradient>(H).eval();
326 const auto E = transformStrain<StrainTags::displacementGradient, StrainTags::greenLagrangian>(H).eval();
328 const auto p = pFunction.evaluate(local, Dune::on(gridElement)).eval();
330 auto [dev, vol] = [&]() {
331 if constexpr (isSameResultType<RT, ResultTypes::PK2StressFull> and PRE::Material::isReduced)
332 return std::make_pair(constitutiveDriver_.deviatoricMaterial().underlying(),
333 constitutiveDriver_.volumetricMaterial().underlying());
335 return std::make_pair(constitutiveDriver_.deviatoricMaterial(), constitutiveDriver_.volumetricMaterial());
341 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>)
342 resultWrapper = PK2Stress;
343 else if constexpr (isSameResultType<RT, ResultTypes::kirchhoffStress> or
344 isSameResultType<RT, ResultTypes::cauchyStress>) {
345 const auto PK2StressMat =
fromVoigt(PK2Stress,
false).eval();
348 resultWrapper =
toVoigt(transformStress<StressTags::PK2, toStress>(PK2StressMat, F),
false);
350 return resultWrapper;
355 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
356 auto& underlying() {
return static_cast<FE&
>(*this); }
357 std::shared_ptr<const Geometry> geo_;
358 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisTypeU>> localBasisU_;
359 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisTypeP>> localBasisP_;
361 size_t numberOfNodes_{0};
375 template <
typename ST>
378 return [&](
const FE::template
KgType<ST>& kgIJ,
const int I,
const int J,
const auto& gp) {
379 const auto geo = underlying().localView().element().geometry();
380 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
381 K.template block<FE::myDim, FE::myDim>(I *
FE::myDim, J *
FE::myDim) += kgIJ * intElement;
394 template <
typename ST>
398 [&](
const auto& C,
const BopType<ST>& bopI,
const BopType<ST>& bopJ,
const int I,
const int J,
const auto& gp) {
399 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
400 K.template block<myDim, myDim>(I *
myDim, J *
myDim) += (bopI.transpose() * C * bopJ) * intElement;
415 template <
typename ST>
419 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
420 force.template segment<myDim>(
myDim * I) += bopI.transpose() * stresses * intElement;
432 template <
typename ST>
435 using namespace Dune::DerivativeDirections;
436 using namespace Dune;
440 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
441 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
442 const auto p = pFunction.evaluate(gpIndex, on(gridElement)).eval();
444 auto e = constitutiveDriver<ST>().storedEnergy(EVoigt, p[0]);
445 energy += e * geo_->integrationElement(gp.position()) * gp.weight();
459 template <
typename ScalarType>
461 typename Traits::template MatrixType<> K,
464 static_assert(Dune::AlwaysFalse<ScalarType>::value,
465 "DisplacementPressure element does not support matrix calculations for autodiff scalars");
467 using namespace Dune::DerivativeDirections;
468 using namespace Dune;
471 auto pFE = underlying().localView().tree().child(Dune::Indices::_1);
472 auto plocalBasis = pFE.finiteElement().localBasis();
473 std::vector<Dune::FieldVector<double, 1>> Np;
476 const auto kMFunction = materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
477 const auto kGFunction = geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
478 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
479 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
480 const auto p = pFunction.evaluate(gpIndex, on(gridElement)).eval();
481 plocalBasis.evaluateFunction(gp.position(), Np);
482 const auto intElement = geo_->integrationElement(gp.position()) * gp.weight();
484 auto [dUdE, S] = constitutiveDriver<ScalarType>().firstDerivatives(EVoigt, p[0]);
485 auto C = constitutiveDriver<ScalarType>().secondDerivative(EVoigt, p[0]);
488 for (
size_t i = 0; i < numberOfNodes_; ++i) {
489 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
490 for (
size_t j = 0; j < numberOfNodes_; ++j) {
491 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
492 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(S), on(gridElement));
493 kMFunction(C, bopI, bopJ, i, j, gp);
494 kGFunction(kgIJ, i, j, gp);
499 for (
size_t i = 0; i < numberOfNodes_; ++i) {
500 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
501 for (
auto j : Dune::range(pFE.size())) {
502 auto jIdx = pFE.localIndex(j);
503 const auto upTerm = (bopI.transpose() * dUdE * Np[j][0]).eval();
504 K.template block<myDim, 1>(i *
myDim, jIdx) += upTerm * intElement;
509 for (
auto i : Dune::range(pFE.size())) {
510 auto iIdx = pFE.localIndex(i);
511 for (
size_t j = 0; j < numberOfNodes_; ++j) {
512 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
513 const auto upTerm = (bopJ.transpose() * dUdE * Np[i][0]).
transpose().eval();
514 K.template block<1, myDim>(iIdx, j *
myDim) += upTerm * intElement;
519 for (
auto i : Dune::range(pFE.size())) {
520 auto iIdx = pFE.localIndex(i);
521 for (
auto j : Dune::range(pFE.size())) {
522 auto jIdx = pFE.localIndex(j);
523 auto ppTerm = Np[i][0] * (1.0 / constitutiveDriver_.kappa()) * Np[j][0];
524 K(iIdx, jIdx) -= ppTerm * intElement;
530 template <
typename ScalarType>
536 static_assert(Dune::AlwaysFalse<ScalarType>::value,
537 "DisplacementPressure element does not support scalar calculations for autodiff scalars");
541 template <
typename ScalarType>
543 typename Traits::template VectorType<ScalarType> force,
546 static_assert(Dune::AlwaysFalse<ScalarType>::value,
547 "DisplacementPressure element does not support vector calculations for autodiff scalars");
549 using namespace Dune::DerivativeDirections;
550 using namespace Dune;
554 auto pFE = underlying().localView().tree().child(Dune::Indices::_1);
555 auto plocalBasis = pFE.finiteElement().localBasis();
556 std::vector<Dune::FieldVector<double, 1>> Np;
558 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
560 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
561 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
562 const auto p = pFunction.evaluate(gpIndex, on(gridElement)).eval();
563 plocalBasis.evaluateFunction(gp.position(), Np);
565 auto [dUdE, S] = constitutiveDriver<ScalarType>().firstDerivatives(EVoigt, p[0]);
566 auto Uhat = constitutiveDriver<ScalarType>().Uhat(EVoigt);
569 for (
size_t i = 0; i < numberOfNodes_; ++i) {
570 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
571 fIntFunction(S, bopI, i, gp);
575 for (
auto i : Dune::range(pFE.size())) {
576 auto iIdx = pFE.localIndex(i);
577 auto pTerm = Np[i][0] * (Uhat - p[0] / constitutiveDriver_.kappa());
578 force(iIdx) += pTerm * geo_->integrationElement(gp.position()) * gp.weight();
590template <Concepts::Material MAT>
599 #error DisplacementPressure depends on dune-localfefunctions, which is not included
Helper for transform between Dune linear algebra types and Eigen.
Helper for the autodiff library.
Free function for creating a tensor product quadrature rule and other Dune::Quadrature rule utils.
Collection of fallback default functions.
Material property functions and conversion utilities.
Header file for types of loads in Ikarus finite element mechanics.
Implementation of transformations for different strain measures.
Implementation of transformations for different stress measures.
Definition of several material related enums.
Definition of helper struct and function to decompose a hyperelastic material model.
Header file for various EAS functions.
Definition of the LinearElastic class for finite element mechanics computations.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
StressTags
A strongly typed enum class representing the type of the computed stresses.
Definition: tags.hh:26
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 displacementPressure(const MAT &mat)
A helper function to create a displacement-pressure pre finite element.
Definition: displacementpressure.hh:591
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 transpose(const Eigen::EigenBase< Derived > &A)
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:39
auto decomposeHyperelasticAndGetMaterialParameters(const MAT &mat)
A helper function to decompose a hyperelastic material model and get the underlying deviatoric materi...
Definition: decomposehyperelastic.hh:47
Definition: utils/dirichletvalues.hh:36
FE class is a base class for all finite elements.
Definition: febase.hh:79
static constexpr int myDim
Definition: febase.hh:95
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
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:278
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
DisplacementPressure class represents a displacement-pressure finite element.
Definition: displacementpressure.hh:133
decltype(std::declval< LocalView >().tree().child(Dune::Indices::_0).child(0).finiteElement().localBasis()) LocalBasisTypeU
Definition: displacementpressure.hh:146
static constexpr auto strainType
Definition: displacementpressure.hh:155
auto calculateStress(const D &dev, const V &vol, const Eigen::Vector< double, strainDim > &strainInVoigt, double p) const
Get the PK2 stress for a given strain (in Voigt notation).
Definition: displacementpressure.hh:300
DisplacementPressure(const Pre &pre)
Constructor for the DisplacementPressure class.
Definition: displacementpressure.hh:185
typename Traits::Geometry Geometry
Definition: displacementpressure.hh:140
static constexpr auto stressType
Definition: displacementpressure.hh:156
decltype(auto) constitutiveDriver() const
Definition: displacementpressure.hh:277
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: displacementpressure.hh:318
typename Pre::Material MaterialType
Definition: displacementpressure.hh:170
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: displacementpressure.hh:542
int order() const
Definition: displacementpressure.hh:272
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: displacementpressure.hh:151
typename Traits::Basis Basis
Definition: displacementpressure.hh:136
auto geometricStiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the geometric part of the stiffness matrix (Kg) for a given inte...
Definition: displacementpressure.hh:376
typename Traits::Element Element
Definition: displacementpressure.hh:142
Impl::DPConstitutiveDriver< DecomposedDevType< MaterialType >, DecomposedVolType< MaterialType > > DPConstitutiveDriverType
Definition: displacementpressure.hh:179
typename Materials::DecomposedMaterialTypes< MaterialType >::DevType DecomposedDevType
Definition: displacementpressure.hh:173
const Geometry & geometry() const
Definition: displacementpressure.hh:270
size_t numberOfNodes() const
Definition: displacementpressure.hh:271
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisTypeU > > & localBasis() const
Definition: displacementpressure.hh:273
typename Traits::GridView GridView
Definition: displacementpressure.hh:141
typename Traits::LocalView LocalView
Definition: displacementpressure.hh:139
void bindImpl()
A helper function to bind the local view to the element.
Definition: displacementpressure.hh:197
auto materialStiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the material part of the stiffness matrix (Ke + Ku) for a given ...
Definition: displacementpressure.hh:395
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: displacementpressure.hh:168
PRE Pre
Definition: displacementpressure.hh:143
typename Materials::DecomposedMaterialTypes< MaterialType >::VolType DecomposedVolType
Definition: displacementpressure.hh:176
static constexpr int strainDim
Definition: displacementpressure.hh:154
typename Traits::FlatBasis FlatBasis
Definition: displacementpressure.hh:137
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Calculate the matrix associated with the given Requirement.
Definition: displacementpressure.hh:460
Eigen::Vector< ST, strainDim > StrainType
Definition: displacementpressure.hh:162
auto strainFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
The strain function for the given Requirement.
Definition: displacementpressure.hh:266
auto displacementFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get the displacement function for the given Requirement.
Definition: displacementpressure.hh:231
auto internalForcesFunction(const Requirement &par, typename Traits::template VectorType< ST > &force, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the internal force vector for a given strain,...
Definition: displacementpressure.hh:416
decltype(std::declval< LocalView >().tree().child(Dune::Indices::_1).finiteElement().localBasis()) LocalBasisTypeP
Definition: displacementpressure.hh:148
static constexpr int myDim
Definition: displacementpressure.hh:153
auto pressureFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get the pressure function for the given Requirement.
Definition: displacementpressure.hh:248
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: displacementpressure.hh:165
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: displacementpressure.hh:531
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisTypeP > > & localBasisP() const
Definition: displacementpressure.hh:274
auto energyFunction(const Requirement &par, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the internal energy at a given integration point and its index.
Definition: displacementpressure.hh:433
A PreFE struct for displacement-pressure elements.
Definition: displacementpressure.hh:113
MAT material
Definition: displacementpressure.hh:115
MAT Material
Definition: displacementpressure.hh:114
std::tuple_element_t< 1, TupleType > VolType
Type of the pure volumetric material model.
Definition: decomposehyperelastic.hh:25
std::tuple_element_t< 0, TupleType > DevType
Type of the pure deviatoric material model.
Definition: decomposehyperelastic.hh:24
Definition: utils/dirichletvalues.hh:38
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:625
Implementation of the Hyperelastic material model.