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>
34template <
typename PreFE,
typename FE,
typename PRE>
35class NonLinearElastic;
41template <Concepts::Material MAT>
47 template <
typename PreFE,
typename FE>
60template <
typename PreFE,
typename FE,
typename PRE>
75 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
77 template <
typename ST>
78 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
85 template <
template <
typename,
int,
int>
class RT>
88 template <
typename ST =
double>
91 template <
typename ST =
double>
92 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
94 template <
typename ST =
double>
95 using KgType = Eigen::Matrix<ST, myDim, myDim>;
109 const auto& localView = underlying().localView();
110 const auto& element = localView.element();
111 auto& firstChild = localView.tree().child(0);
112 const auto& fe = firstChild.finiteElement();
113 geo_ = std::make_shared<const Geometry>(element.geometry());
114 numberOfNodes_ = fe.size();
115 order_ = 2 * (fe.localBasis().order());
116 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
117 if constexpr (
requires { element.impl().getQuadratureRule(order_); })
118 if (element.impl().isTrimmed())
119 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
121 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
122 Dune::bindDerivatives(0, 1));
124 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
136 template <
typename ScalarType =
double>
139 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
140 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
152 template <
typename ScalarType =
double>
165 template <
typename ScalarType,
int strainDim>
167 return material<ScalarType>().template storedEnergy<strainType>(strain);
179 template <
typename ScalarType,
int strainDim,
bool voigt = true>
180 auto stress(
const Eigen::Vector<ScalarType, strainDim>& strain)
const {
181 return material<ScalarType>().template stresses<strainType, voigt>(strain);
193 template <
typename ScalarType,
int strainDim,
bool voigt = true>
195 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
200 [[nodiscard]]
int order()
const {
return order_; }
201 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>&
localBasis()
const {
return localBasis_; }
203 template <
typename ScalarType =
double>
206 return mat_.template rebind<ScalarType>();
217 template <
template <
typename,
int,
int>
class RT>
218 requires(supportsResultType<RT>())
220 return [&](
const Eigen::Vector<double, strainDim>& strainInVoigt) {
221 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>) {
222 decltype(
auto) mat = [&]() {
223 if constexpr (isSameResultType<RT, ResultTypes::PK2StressFull> and Material::isReduced)
224 return mat_.underlying();
228 return RTWrapperType<RT>{mat.template stresses<strainType>(enlargeIfReduced<Material>(strainInVoigt))};
242 template <
template <
typename,
int,
int>
class RT>
243 requires(supportsResultType<RT>())
245 Dune::PriorityTag<1>)
const {
246 using namespace Dune::DerivativeDirections;
250 if constexpr (isSameResultType<RT, ResultTypes::PK2Stress> or isSameResultType<RT, ResultTypes::PK2StressFull>) {
252 const auto rFunction = resultFunction<RT>();
253 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(gridElement));
254 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
262 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
263 auto& underlying() {
return static_cast<FE&
>(*this); }
264 std::shared_ptr<const Geometry> geo_;
265 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
267 size_t numberOfNodes_{0};
281 template <
typename ST>
284 return [&](
const FE::template
KgType<ST>& kgIJ,
const int I,
const int J,
const auto& gp) {
285 const auto geo = underlying().localView().element().geometry();
286 const double intElement = geo.integrationElement(gp.position()) * gp.weight();
287 K.template block<FE::myDim, FE::myDim>(I *
FE::myDim, J *
FE::myDim) += kgIJ * intElement;
300 template <
typename ST>
306 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
307 K.template block<myDim, myDim>(I *
myDim, J *
myDim) += (bopI.transpose() * C * bopJ) * intElement;
322 template <
typename ST>
326 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
327 force.template segment<myDim>(
myDim * I) += bopI.transpose() * stresses * intElement;
339 template <
typename ST>
342 using namespace Dune::DerivativeDirections;
343 using namespace Dune;
346 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
347 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
348 energy +=
internalEnergy(EVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
362 template <
typename ScalarType>
364 typename Traits::template MatrixType<> K,
368 using namespace Dune::DerivativeDirections;
369 using namespace Dune;
372 const auto kMFunction = materialStiffnessMatrixFunction<ScalarType>(par, K, dx);
373 const auto kGFunction = geometricStiffnessMatrixFunction<ScalarType>(par, K, dx);
374 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
375 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
376 const auto stresses =
stress(EVoigt);
377 for (
size_t i = 0; i < numberOfNodes_; ++i) {
378 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
379 for (
size_t j = 0; j < numberOfNodes_; ++j) {
380 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
381 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(gridElement));
382 kMFunction(EVoigt, bopI, bopJ, i, j, gp);
383 kGFunction(kgIJ, i, j, gp);
389 template <
typename ScalarType>
393 return ScalarType{0.0};
397 template <
typename ScalarType>
399 typename Traits::template VectorType<ScalarType> force,
403 using namespace Dune::DerivativeDirections;
404 using namespace Dune;
406 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
409 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
410 const auto EVoigt = (eps.evaluate(gpIndex, on(gridElement))).eval();
411 const auto stresses =
stress(EVoigt);
412 for (
size_t i = 0; i < numberOfNodes_; ++i) {
413 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
414 fIntFunction(stresses, bopI, i, gp);
426template <Concepts::Material MAT>
436 #error NonLinearElastic depends on dune-localfefunctions, which is not included
Collection of fallback default functions.
Helper for the autodiff library.
Helper for transform between Dune linear algebra types and Eigen.
Definition of several material related enums.
Header file for types of loads in Ikarus finite element mechanics.
Header file for various EAS functions.
Definition of the LinearElastic class for finite element mechanics computations.
Material property functions and conversion utilities.
Contains the FE class, which is used as a base class for all finite elements. It provides information...
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
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:64
auto nonLinearElastic(const MAT &mat)
A helper function to create a non-linear elastic pre finite element.
Definition: nonlinearelastic.hh:427
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:49
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:38
Definition: utils/dirichletvalues.hh:30
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:223
const SolutionVectorType & globalSolution() const
Get the global solution vector.
Definition: ferequirements.hh:313
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:164
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:277
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
NonLinearElastic class represents a non-linear elastic finite element.
Definition: nonlinearelastic.hh:62
auto strainFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
The strain function for the given Requirement.
Definition: nonlinearelastic.hh:153
typename Traits::Basis Basis
Definition: nonlinearelastic.hh:65
typename Traits::Geometry Geometry
Definition: nonlinearelastic.hh:69
typename Traits::FlatBasis FlatBasis
Definition: nonlinearelastic.hh:66
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: nonlinearelastic.hh:92
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: nonlinearelastic.hh:340
PRE::Material Material
Definition: nonlinearelastic.hh:72
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > & localBasis() const
Definition: nonlinearelastic.hh:201
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: nonlinearelastic.hh:78
typename Traits::GridView GridView
Definition: nonlinearelastic.hh:70
void bindImpl()
A helper function to bind the local view to the element.
Definition: nonlinearelastic.hh:108
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: nonlinearelastic.hh:398
const Geometry & geometry() const
Definition: nonlinearelastic.hh:198
static constexpr int myDim
Definition: nonlinearelastic.hh:80
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: nonlinearelastic.hh:390
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: nonlinearelastic.hh:244
static constexpr auto stressType
Definition: nonlinearelastic.hh:83
decltype(auto) material() const
Definition: nonlinearelastic.hh:204
NonLinearElastic(const Pre &pre)
Constructor for the NonLinearElastic class.
Definition: nonlinearelastic.hh:101
typename Traits::Element Element
Definition: nonlinearelastic.hh:71
PRE Pre
Definition: nonlinearelastic.hh:73
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: nonlinearelastic.hh:95
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: nonlinearelastic.hh:282
auto resultFunction() const
Get a lambda function that evaluates the requested result type for a given strain (in Voigt notation)...
Definition: nonlinearelastic.hh:219
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the material tangent for the given strain for the given Requirement.
Definition: nonlinearelastic.hh:194
static constexpr int strainDim
Definition: nonlinearelastic.hh:81
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: nonlinearelastic.hh:301
auto internalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: nonlinearelastic.hh:166
static constexpr auto strainType
Definition: nonlinearelastic.hh:82
Eigen::Vector< ST, strainDim > StrainType
Definition: nonlinearelastic.hh:89
typename Traits::LocalView LocalView
Definition: nonlinearelastic.hh:68
int order() const
Definition: nonlinearelastic.hh:200
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: nonlinearelastic.hh:75
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: nonlinearelastic.hh:323
auto stress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: nonlinearelastic.hh:180
auto displacementFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get the displacement function for the given Requirement.
Definition: nonlinearelastic.hh:137
size_t numberOfNodes() const
Definition: nonlinearelastic.hh:199
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: nonlinearelastic.hh:363
A PreFE struct for non-linear elastic elements.
Definition: nonlinearelastic.hh:43
MAT Material
Definition: nonlinearelastic.hh:44
MAT material
Definition: nonlinearelastic.hh:45
static constexpr bool isMixed()
Definition: mixin.hh:110
Definition: utils/dirichletvalues.hh:32
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:622