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>
29template <
typename PreFE,
typename FE,
typename PRE>
36template <Concepts::GeometricallyLinearMaterial MAT>
42 template <
typename PreFE,
typename FE>
54template <
typename PreFE,
typename FE,
typename PRE>
71 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
78 : mat_{pre.material} {}
85 const auto& localView = underlying().localView();
86 const auto& element = localView.element();
87 auto& firstChild = localView.tree().child(0);
88 const auto& fe = firstChild.finiteElement();
89 geo_ = std::make_shared<const Geometry>(element.geometry());
90 numberOfNodes_ = fe.size();
91 order_ = 2 * (fe.localBasis().order());
92 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
93 if constexpr (
requires { element.impl().getQuadratureRule(order_); })
94 if (element.impl().isTrimmed())
95 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
97 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
98 Dune::bindDerivatives(0, 1));
100 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
112 template <
typename ScalarType =
double>
115 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
117 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
118 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
130 template <
class ScalarType =
double>
133 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
144 return mat_.template tangentModuli<StrainTags::linear, true>(Eigen::Vector<double, 6>::Zero());
159 [[nodiscard]]
int order()
const {
return order_; }
172 template <
template <
typename,
int,
int>
class RT>
173 requires(supportsResultType<RT>())
175 Dune::PriorityTag<1>)
const {
177 if constexpr (isSameResultType<RT, ResultTypes::linearStress>) {
180 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
182 return RTWrapper{(C * epsVoigt).eval()};
188 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
189 auto& underlying() {
return static_cast<FE&
>(*this); }
191 std::shared_ptr<const Geometry> geo_;
192 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
194 size_t numberOfNodes_{0};
198 template <
typename ScalarType>
201 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
203 using namespace Dune::DerivativeDirections;
204 using namespace Dune;
207 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
208 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
209 for (
size_t i = 0; i < numberOfNodes_; ++i) {
210 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
211 for (
size_t j = 0; j < numberOfNodes_; ++j) {
212 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
213 K.template block<myDim, myDim>(i *
myDim, j *
myDim) += bopI.transpose() * C * bopJ * intElement;
219 template <
typename ScalarType>
221 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx =
222 std::nullopt)
const -> ScalarType {
225 const auto& lambda = par.parameter();
226 using namespace Dune::DerivativeDirections;
227 using namespace Dune;
231 ScalarType energy = 0.0;
232 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
233 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
234 energy += (0.5 * EVoigt.dot(C * EVoigt)) * geo_->integrationElement(gp.position()) * gp.weight();
239 template <
typename ScalarType>
242 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
244 using namespace Dune::DerivativeDirections;
245 using namespace Dune;
250 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
251 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
252 auto stresses = (C * eps.evaluate(gpIndex, on(gridElement))).eval();
253 for (
size_t i = 0; i < numberOfNodes_; ++i) {
254 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
255 force.template segment<myDim>(
myDim * i) += bopI.transpose() * stresses * intElement;
267template <
typename MAT>
276 #error LinearElastic depends on dune-localfefunctions, which is not included
Definitions of ResultTypes used for finite element results.
Material property functions and conversion utilities.
Definition of the LinearElastic class for finite element mechanics computations.
Definition: assemblermanipulatorbuildingblocks.hh:22
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
auto linearElastic(const MAT &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:268
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:252
SolutionVectorReturnType globalSolution()
Get the global solution vector.
Definition: ferequirements.hh:308
Container that is used for FE Results. It gives access to the stored value, but can also be used to a...
Definition: feresulttypes.hh:159
Base class for element definitions that provides common functionality for ResultTypes.
Definition: feresulttypes.hh:272
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:56
void bindImpl()
A helper function to bind the local view to the element.
Definition: linearelastic.hh:84
static constexpr int myDim
Definition: linearelastic.hh:70
auto materialTangent() const
Gets the material tangent matrix for the linear elastic material.
Definition: linearelastic.hh:142
typename Traits::Element Element
Definition: linearelastic.hh:66
auto materialTangentFunction(const Requirement &par) const
Gets the material tangent function for the given Requirement.
Definition: linearelastic.hh:153
PRE Pre
Definition: linearelastic.hh:68
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const -> ScalarType
Definition: linearelastic.hh:220
size_t numberOfNodes() const
Definition: linearelastic.hh:158
int order() const
Definition: linearelastic.hh:159
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:174
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ScalarType > force, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Definition: linearelastic.hh:240
typename Traits::GridView GridView
Definition: linearelastic.hh:65
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:59
typename Traits::LocalView LocalView
Definition: linearelastic.hh:63
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:60
void calculateMatrixImpl(const Requirement &par, const MatrixAffordance &affordance, typename Traits::template MatrixType<> K, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Definition: linearelastic.hh:199
auto displacementFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Gets the displacement function for the given Requirement and optional displacement vector.
Definition: linearelastic.hh:113
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:77
PRE::Material Material
Definition: linearelastic.hh:67
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:71
typename Traits::Geometry Geometry
Definition: linearelastic.hh:64
const Geometry & geometry() const
Definition: linearelastic.hh:157
auto strainFunction(const Requirement &par, const std::optional< std::reference_wrapper< const Eigen::VectorX< ScalarType > > > &dx=std::nullopt) const
Gets the strain function for the given Requirement and optional di splacement vector.
Definition: linearelastic.hh:131
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:38
MAT Material
Definition: linearelastic.hh:39
MAT material
Definition: linearelastic.hh:40
Definition: utils/dirichletvalues.hh:32
Header file for material models in Ikarus finite element mechanics.