12#if HAVE_DUNE_LOCALFEFUNCTIONS
15 #include <type_traits>
17 #include <dune/geometry/quadraturerules.hh>
18 #include <dune/localfefunctions/expressions/linearStrainsExpr.hh>
19 #include <dune/localfefunctions/impl/standardLocalFunction.hh>
31template <
typename PreFE,
typename FE,
typename PRE>
38template <Concepts::GeometricallyLinearMaterial MAT>
44 template <
typename PreFE,
typename FE>
56template <
typename PreFE,
typename FE,
typename PRE>
71 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
73 template <
typename ST>
74 using VectorXOptRef = std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>;
80 static constexpr bool hasEAS = FE::template hasEAS<EAS::LinearStrain>;
82 template <
template <
typename,
int,
int>
class RT>
85 template <
typename ST>
88 template <
typename ST>
89 using BopType = Eigen::Matrix<ST, strainDim, myDim>;
91 template <
typename ST>
92 using KgType = Eigen::Matrix<ST, myDim, myDim>;
106 const auto& localView = underlying().localView();
107 const auto& element = localView.element();
108 auto& firstChild = localView.tree().child(0);
109 const auto& fe = firstChild.finiteElement();
110 geo_ = std::make_shared<const Geometry>(element.geometry());
111 numberOfNodes_ = fe.size();
112 order_ = 2 * (fe.localBasis().order());
113 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
114 if constexpr (
requires { element.impl().getQuadratureRule(order_); })
115 if (element.impl().isTrimmed())
116 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
118 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
119 Dune::bindDerivatives(0, 1));
121 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
133 template <
typename ScalarType =
double>
136 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
137 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
149 template <
class ScalarType =
double>
162 template <
typename ScalarType,
int strainDim>
164 const auto C = materialTangent<ScalarType, strainDim>(strain);
165 return (0.5 * strain.dot(C * strain));
177 template <
typename ScalarType,
int strainDim,
bool voigt = true>
178 auto stress(
const Eigen::Vector<ScalarType, strainDim>& strain)
const {
179 const auto C = materialTangent<ScalarType, strainDim, voigt>(strain);
180 return (C * strain).eval();
192 template <
typename ScalarType,
int strainDim,
bool voigt = true>
194 const Eigen::Vector<ScalarType, strainDim>& strain = Eigen::Vector<ScalarType, strainDim>::Zero())
const {
196 return material<ScalarType>().template tangentModuli<strainType, voigt>(strain);
201 [[nodiscard]]
int order()
const {
return order_; }
202 const Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>&
localBasis()
const {
return localBasis_; }
204 template <
typename ScalarType =
double>
207 return mat_.template rebind<ScalarType>();
218 template <
template <
typename,
int,
int>
class RT>
219 requires(LinearElastic::template supportsResultType<RT>())
221 return [&](
const Eigen::Vector<double, strainDim>& strainInVoigt) {
222 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or
223 isSameResultType<RT, ResultTypes::linearStressFull>) {
224 decltype(
auto) mat = [&]() {
225 if constexpr (isSameResultType<RT, ResultTypes::linearStressFull> and
requires { mat_.underlying(); })
226 return mat_.underlying();
230 return RTWrapperType<RT>{mat.template stresses<StrainTags::linear>(enlargeIfReduced<Material>(strainInVoigt))};
243 template <
template <
typename,
int,
int>
class RT>
244 requires(LinearElastic::template supportsResultType<RT>())
246 Dune::PriorityTag<1>)
const {
249 if constexpr (isSameResultType<RT, ResultTypes::linearStress> or
250 isSameResultType<RT, ResultTypes::linearStressFull>) {
251 const auto rFunction = resultFunction<RT>();
253 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
254 return rFunction(epsVoigt);
260 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
261 auto& underlying() {
return static_cast<FE&
>(*this); }
263 std::shared_ptr<const Geometry> geo_;
264 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
266 size_t numberOfNodes_{0};
280 template <
typename ST>
284 const int I,
const int J,
const auto& gp) {
286 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
287 K.template block<myDim, myDim>(I *
myDim, J *
myDim) += (bopI.transpose() * C * bopJ) * intElement;
302 template <
typename ST>
306 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
307 auto stresses =
stress(strain);
308 force.template segment<myDim>(
myDim * I) += bopI.transpose() * stresses * intElement;
320 template <
typename ScalarType>
322 return [&]() -> ScalarType {
323 using namespace Dune::DerivativeDirections;
324 using namespace Dune;
325 ScalarType energy = 0.0;
327 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
328 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
329 energy +=
internalEnergy(epsVoigt) * geo_->integrationElement(gp.position()) * gp.weight();
336 template <
typename ScalarType>
338 typename Traits::template MatrixType<> K,
342 using namespace Dune::DerivativeDirections;
343 using namespace Dune;
345 const auto kFunction = stiffnessMatrixFunction<ScalarType>(par, K, dx);
348 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
349 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
350 for (
size_t i = 0; i < numberOfNodes_; ++i) {
351 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
352 for (
size_t j = 0; j < numberOfNodes_; ++j) {
353 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
354 kFunction(epsVoigt, bopI, bopJ, kgIJ, i, j, gp);
360 template <
typename ScalarType>
364 return ScalarType{0.0};
368 template <
typename ScalarType>
370 typename Traits::template VectorType<ScalarType> force,
374 using namespace Dune::DerivativeDirections;
375 using namespace Dune;
377 const auto fIntFunction = internalForcesFunction<ScalarType>(par, force, dx);
380 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
381 const auto epsVoigt = eps.evaluate(gpIndex, on(gridElement));
382 for (
size_t i = 0; i < numberOfNodes_; ++i) {
383 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
384 fIntFunction(epsVoigt, bopI, i, gp);
396template <
typename MAT>
405 #error LinearElastic depends on dune-localfefunctions, which is not included
Definition of the LinearElastic class for finite element mechanics computations.
Definitions of ResultTypes used for finite element results.
Material property functions and conversion utilities.
Definition of several material related enums.
Header file for various EAS functions.
Definition: assemblermanipulatorbuildingblocks.hh:22
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:64
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
auto linearElastic(const MAT &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:397
Definition: utils/dirichletvalues.hh:30
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: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
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:58
void bindImpl()
A helper function to bind the local view to the element.
Definition: linearelastic.hh:105
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: linearelastic.hh:303
auto energyFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Get a lambda function that evaluates the internal energy at a given integration point and its index.
Definition: linearelastic.hh:321
static constexpr int myDim
Definition: linearelastic.hh:76
Eigen::Vector< ST, strainDim > StrainType
Definition: linearelastic.hh:86
Eigen::Matrix< ST, strainDim, myDim > BopType
Definition: linearelastic.hh:89
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: linearelastic.hh:337
Eigen::Matrix< ST, myDim, myDim > KgType
Definition: linearelastic.hh:92
typename Traits::Element Element
Definition: linearelastic.hh:67
auto resultFunction() const
Get a lambda function that evaluates the requested result type for a given strain (in Voigt notation)...
Definition: linearelastic.hh:220
static constexpr auto strainType
Definition: linearelastic.hh:78
auto stiffnessMatrixFunction(const Requirement &par, typename Traits::template MatrixType< ST > &K, const VectorXOptRef< ST > &dx=std::nullopt) const
Get a lambda function that evaluates the stiffness matrix for a given strain, integration point and i...
Definition: linearelastic.hh:281
PRE Pre
Definition: linearelastic.hh:69
static constexpr auto stressType
Definition: linearelastic.hh:79
static constexpr bool hasEAS
Definition: linearelastic.hh:80
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const VectorXOptRef< ScalarType > &dx=std::nullopt) const -> ScalarType
Definition: linearelastic.hh:361
size_t numberOfNodes() const
Definition: linearelastic.hh:200
int order() const
Definition: linearelastic.hh:201
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain=Eigen::Vector< ScalarType, strainDim >::Zero()) const
Get the linear elastic material tangent for the given strain for the given Requirement.
Definition: linearelastic.hh:193
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:245
typename Traits::GridView GridView
Definition: linearelastic.hh:66
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:61
typename Traits::LocalView LocalView
Definition: linearelastic.hh:64
std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > VectorXOptRef
Definition: linearelastic.hh:74
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:62
auto displacementFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Gets the displacement function for the given Requirement and optional displacement vector.
Definition: linearelastic.hh:134
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:98
PRE::Material Material
Definition: linearelastic.hh:68
static constexpr int strainDim
Definition: linearelastic.hh:77
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:71
decltype(auto) material() const
Definition: linearelastic.hh:205
typename Traits::Geometry Geometry
Definition: linearelastic.hh:65
auto internalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: linearelastic.hh:163
const Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > & localBasis() const
Definition: linearelastic.hh:202
auto strainFunction(const Requirement &par, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Gets the strain function for the given Requirement and optional di splacement vector.
Definition: linearelastic.hh:150
const Geometry & geometry() const
Definition: linearelastic.hh:199
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const VectorXOptRef< ScalarType > &dx=std::nullopt) const
Definition: linearelastic.hh:369
auto stress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: linearelastic.hh:178
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:40
MAT Material
Definition: linearelastic.hh:41
MAT material
Definition: linearelastic.hh:42
Definition: utils/dirichletvalues.hh:32
Concept to check if the underlying scalar type is a dual type.
Definition: utils/concepts.hh:615
Header file for material models in Ikarus finite element mechanics.