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>
28template <
typename PreFE,
typename FE>
38 template <
typename PreFE,
typename FE>
50template <
typename PreFE,
typename FE>
66 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
73 : mat_{pre.material} {}
80 const auto& localView = underlying().localView();
81 const auto& element = localView.element();
82 auto& firstChild = localView.tree().child(0);
83 const auto& fe = firstChild.finiteElement();
84 geo_ = std::make_shared<const Geometry>(element.geometry());
85 numberOfNodes_ = fe.size();
86 order_ = 2 * (fe.localBasis().order());
87 localBasis_ = Dune::CachedLocalBasis(fe.localBasis());
88 if constexpr (
requires { element.impl().getQuadratureRule(order_); })
89 if (element.impl().isTrimmed())
90 localBasis_.bind(element.impl().getQuadratureRule(order_), Dune::bindDerivatives(0, 1));
92 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_),
93 Dune::bindDerivatives(0, 1));
95 localBasis_.bind(Dune::QuadratureRules<double, myDim>::rule(element.type(), order_), Dune::bindDerivatives(0, 1));
107 template <
typename ScalarType =
double>
110 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
112 auto disp = Ikarus::FEHelper::localSolutionBlockVector<Traits>(d, underlying().localView(), dx);
113 Dune::StandardLocalFunction uFunction(localBasis_, disp, geo_);
124 template <
class ScalarType =
double>
127 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
137 if constexpr (
myDim == 2)
139 else if constexpr (
myDim == 3)
155 [[nodiscard]]
int order()
const {
return order_; }
159 template <
template <
typename,
int,
int>
class RT>
161 return isSameResultType<RT, ResultTypes::linearStress>;
175 template <
template <
typename,
int,
int>
class RT>
176 requires(canProvideResultType<RT>())
178 Dune::PriorityTag<1>)
const {
180 if constexpr (isSameResultType<RT, ResultTypes::linearStress>) {
183 auto epsVoigt = eps.evaluate(local, Dune::on(Dune::DerivativeDirections::gridElement));
185 return RTWrapper{(C * epsVoigt).eval()};
191 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
192 auto& underlying() {
return static_cast<FE&
>(*this); }
194 std::shared_ptr<const Geometry> geo_;
195 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>> localBasis_;
196 YoungsModulusAndPoissonsRatio mat_;
197 size_t numberOfNodes_{0};
201 template <
typename ScalarType>
204 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
206 using namespace Dune::DerivativeDirections;
207 using namespace Dune;
210 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
211 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
212 for (
size_t i = 0; i < numberOfNodes_; ++i) {
213 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
214 for (
size_t j = 0; j < numberOfNodes_; ++j) {
215 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(gridElement));
216 K.template block<myDim, myDim>(i *
myDim, j *
myDim) += bopI.transpose() * C * bopJ * intElement;
222 template <
typename ScalarType>
224 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx =
225 std::nullopt)
const -> ScalarType {
228 const auto& lambda = par.parameter();
229 using namespace Dune::DerivativeDirections;
230 using namespace Dune;
234 ScalarType energy = 0.0;
235 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
236 const auto EVoigt = eps.evaluate(gpIndex, on(gridElement));
237 energy += (0.5 * EVoigt.dot(C * EVoigt)) * geo_->integrationElement(gp.position()) * gp.weight();
242 template <
typename ScalarType>
245 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ScalarType>>>& dx = std::nullopt)
const {
247 using namespace Dune::DerivativeDirections;
248 using namespace Dune;
253 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
254 const double intElement = geo_->integrationElement(gp.position()) * gp.weight();
255 auto stresses = (C * eps.evaluate(gpIndex, on(gridElement))).eval();
256 for (
size_t i = 0; i < numberOfNodes_; ++i) {
257 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(gridElement));
258 force.template segment<myDim>(
myDim * i) += bopI.transpose() * stresses * intElement;
277 #error LinearElastic depends on dune-localfefunctions, which is not included
Definition of the LinearElastic class for finite element mechanics computations.
Header file for material models in Ikarus finite element mechanics.
Material property functions and conversion utilities.
Definition: dirichletbcenforcement.hh:6
MatrixAffordance
A strongly typed enum class representing the matrix affordance.
Definition: ferequirements.hh:63
Eigen::Matrix3d planeStressLinearElasticMaterialTangent(double E, double nu)
Computes the plane stress linear elastic material tangent matrix.
Definition: physicshelper.hh:23
VectorAffordance
A strongly typed enum class representing the vector affordance.
Definition: ferequirements.hh:48
Eigen::Matrix< double, 6, 6 > linearElasticMaterialTangent3D(double E, double nu)
Computes the 3D linear elastic material tangent matrix.
Definition: physicshelper.hh:40
auto linearElastic(const YoungsModulusAndPoissonsRatio &mat)
A helper function to create a linear elastic pre finite element.
Definition: linearelastic.hh:269
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
Definition: utils/dirichletvalues.hh:28
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:155
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:52
typename Traits::FlatBasis FlatBasis
Definition: linearelastic.hh:56
typename Traits::Geometry Geometry
Definition: linearelastic.hh:60
typename Traits::Element Element
Definition: linearelastic.hh:62
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:223
auto materialTangent() const
Gets the material tangent matrix for the linear elastic material.
Definition: linearelastic.hh:136
auto materialTangentFunction(const Requirement &par) const
Gets the material tangent function for the given Requirement.
Definition: linearelastic.hh:149
int order() const
Definition: linearelastic.hh:155
static consteval bool canProvideResultType()
Definition: linearelastic.hh:160
typename Traits::LocalView LocalView
Definition: linearelastic.hh:59
typename Traits::GridView GridView
Definition: linearelastic.hh:61
LinearElastic(const Pre &pre)
Constructor for the LinearElastic class.
Definition: linearelastic.hh:72
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:243
typename Traits::BasisHandler BasisHandler
Definition: linearelastic.hh:55
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 displacement vector.
Definition: linearelastic.hh:125
const Geometry & geometry() const
Definition: linearelastic.hh:153
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: linearelastic.hh:66
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:202
std::tuple< decltype(makeRT< ResultTypes::linearStress >())> SupportedResultTypes
Definition: linearelastic.hh:157
void bindImpl()
A helper function to bind the local view to the element.
Definition: linearelastic.hh:79
static constexpr int myDim
Definition: linearelastic.hh:65
size_t numberOfNodes() const
Definition: linearelastic.hh:154
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:177
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:108
A PreFE struct for linear elastic elements.
Definition: linearelastic.hh:35
YoungsModulusAndPoissonsRatio material
Definition: linearelastic.hh:36
Structure representing Young's modulus and shear modulus.
Definition: physicshelper.hh:53
double emodul
Definition: physicshelper.hh:54
double nu
Definition: physicshelper.hh:55
Definition: utils/dirichletvalues.hh:30