version 0.4.1
traction.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#pragma once
5
6#include <dune/fufem/boundarypatch.hh>
7#include <dune/localfefunctions/derivativetransformators.hh>
8#include <dune/localfefunctions/meta.hh>
9
11
12namespace Ikarus {
13
14template <typename PreFE, typename FE>
15class Traction;
16
21template <typename GV>
23{
24 using GridView = GV;
25 static constexpr int worldDim = GridView::dimensionworld;
26 const BoundaryPatch<GridView>* neumannBoundary;
27
28 std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)> load;
29 using BoundaryPatchType = BoundaryPatch<GridView>;
30
31 template <typename PreFE, typename FE>
33};
34
43template <typename PreFE, typename FE>
45{
46public:
50 using LocalView = typename Traits::LocalView;
51 using GridView = typename Traits::GridView;
52 static constexpr int myDim = Traits::mydim;
53 static constexpr int worldDim = Traits::worlddim;
55
61 explicit Traction(const Pre& pre)
62 : neumannBoundaryLoad_{pre.load},
63 neumannBoundary_{pre.neumannBoundary} {}
64
65protected:
66 template <template <typename, int, int> class RT>
67 requires Dune::AlwaysFalse<RT<double, 1, 1>>::value
69 Dune::PriorityTag<0>) const {}
70
71 template <typename ST>
73 const Requirement& par, ScalarAffordance affordance,
74 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const -> ST {
75 if (not neumannBoundary_ and not neumannBoundaryLoad_)
76 return 0.0;
77 ST energy = 0.0;
78 const auto uFunction = underlying().displacementFunction(par, dx);
79 const auto& lambda = par.parameter();
80 auto& element = underlying().localView().element();
81
82 for (auto&& intersection : intersections(neumannBoundary_->gridView(), element)) {
83 if (not neumannBoundary_->contains(intersection))
84 continue;
85
86 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), underlying().order());
87
88 for (const auto& curQuad : quadLine) {
89 // Local position of the quadrature point
90 const Dune::FieldVector<double, myDim>& quadPos = intersection.geometryInInside().global(curQuad.position());
91
92 const double integrationElement = intersection.geometry().integrationElement(curQuad.position());
93
94 // The value of the local function
95 const auto uVal = uFunction.evaluate(quadPos);
96
97 // Value of the Neumann data at the current position
98 auto neumannValue = neumannBoundaryLoad_(intersection.geometry().global(curQuad.position()), lambda);
99
100 energy -= neumannValue.dot(uVal) * curQuad.weight() * integrationElement;
101 }
102 }
103 return energy;
104 }
105
106 template <typename ST>
108 const Requirement& par, VectorAffordance affordance, typename Traits::template VectorType<ST> force,
109 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& dx = std::nullopt) const {
110 if (not neumannBoundary_ and not neumannBoundaryLoad_)
111 return;
112 using namespace Dune::DerivativeDirections;
113 using namespace Dune;
114 const auto uFunction = underlying().displacementFunction(par, dx);
115 const auto& lambda = par.parameter();
116 auto& element = underlying().localView().element();
117
118 for (auto&& intersection : intersections(neumannBoundary_->gridView(), element)) {
119 if (not neumannBoundary_->contains(intersection))
120 continue;
121
123 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), underlying().order());
124
125 for (const auto& curQuad : quadLine) {
126 const Dune::FieldVector<double, myDim>& quadPos = intersection.geometryInInside().global(curQuad.position());
127 const double intElement = intersection.geometry().integrationElement(curQuad.position());
128
130 for (size_t i = 0; i < underlying().numberOfNodes(); ++i) {
131 const auto udCi = uFunction.evaluateDerivative(quadPos, wrt(coeff(i)));
132
134 auto neumannValue = neumannBoundaryLoad_(intersection.geometry().global(curQuad.position()), lambda);
135 force.template segment<worldDim>(worldDim * i) -= udCi * neumannValue * curQuad.weight() * intElement;
136 }
137 }
138 }
139 }
140
141 template <typename ST>
143 const Requirement&, MatrixAffordance, typename Traits::template MatrixType<>,
144 const std::optional<std::reference_wrapper<const Eigen::VectorX<ST>>>& = std::nullopt) const {}
145
146private:
147 std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)>
148 neumannBoundaryLoad_;
149 const BoundaryPatch<GridView>* neumannBoundary_;
150
151 //> CRTP
152 const auto& underlying() const { return static_cast<const FE&>(*this); }
153 auto& underlying() { return static_cast<FE&>(*this); }
154};
155
164template <typename GV, typename F>
165auto neumannBoundaryLoad(const BoundaryPatch<GV>* patch, F&& load) {
166 NeumannBoundaryLoadPre<GV> pre(patch, std::forward<F>(load));
167
168 return pre;
169}
170
171} // namespace Ikarus
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
auto neumannBoundaryLoad(const BoundaryPatch< GV > *patch, F &&load)
A helper function to create a Neumann boundary load skill.
Definition: traction.hh:165
ScalarAffordance
A strongly typed enum class representing the scalar affordance.
Definition: ferequirements.hh:37
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
PMHelper::ConstReturnType parameter() const
Get the parameter value.
Definition: ferequirements.hh:325
Traits for handling finite elements.
Definition: fetraits.hh:25
typename Basis::LocalView LocalView
Type of the local view.
Definition: fetraits.hh:42
typename Basis::GridView GridView
Type of the grid view.
Definition: fetraits.hh:45
static constexpr int worlddim
Dimension of the world space.
Definition: fetraits.hh:60
static constexpr int mydim
Dimension of the geometry.
Definition: fetraits.hh:63
Traction class represents distributed traction load that can be applied.
Definition: traction.hh:45
typename Traits::LocalView LocalView
Definition: traction.hh:50
typename Traits::GridView GridView
Definition: traction.hh:51
static constexpr int myDim
Definition: traction.hh:52
static constexpr int worldDim
Definition: traction.hh:53
void calculateVectorImpl(const Requirement &par, VectorAffordance affordance, typename Traits::template VectorType< ST > force, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const
Definition: traction.hh:107
void calculateMatrixImpl(const Requirement &, MatrixAffordance, typename Traits::template MatrixType<>, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &=std::nullopt) const
Definition: traction.hh:142
auto calculateScalarImpl(const Requirement &par, ScalarAffordance affordance, const std::optional< std::reference_wrapper< const Eigen::VectorX< ST > > > &dx=std::nullopt) const -> ST
Definition: traction.hh:72
Traction(const Pre &pre)
Constructor for the traction class.
Definition: traction.hh:61
::value auto calculateAtImpl(const Requirement &req, const Dune::FieldVector< double, Traits::mydim > &local, Dune::PriorityTag< 0 >) const
Definition: traction.hh:68
A PreFE struct for Neumann boundary load skill.
Definition: traction.hh:23
static constexpr int worldDim
Definition: traction.hh:25
GV GridView
Definition: traction.hh:24
const BoundaryPatch< GridView > * neumannBoundary
Definition: traction.hh:26
BoundaryPatch< GridView > BoundaryPatchType
Definition: traction.hh:29
std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> load
Definition: traction.hh:28
Definition: utils/dirichletvalues.hh:32