6#include <dune/fufem/boundarypatch.hh>
7#include <dune/localfefunctions/derivativetransformators.hh>
8#include <dune/localfefunctions/meta.hh>
14template <
typename PreFE,
typename FE>
25 static constexpr int worldDim = GridView::dimensionworld;
31 template <
typename PreFE,
typename FE>
43template <
typename PreFE,
typename FE>
62 : neumannBoundaryLoad_{pre.load},
63 neumannBoundary_{pre.neumannBoundary} {}
66 template <
template <
typename,
int,
int>
class RT>
67 requires Dune::AlwaysFalse<RT<double, 1, 1>>::value
69 Dune::PriorityTag<0>)
const {}
71 template <
typename ST>
74 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const -> ST {
75 if (not neumannBoundary_ and not neumannBoundaryLoad_)
78 const auto uFunction = underlying().displacementFunction(par, dx);
79 const auto& lambda = par.parameter();
80 auto& element = underlying().localView().element();
82 for (
auto&& intersection : intersections(neumannBoundary_->gridView(), element)) {
83 if (not neumannBoundary_->contains(intersection))
86 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), underlying().order());
88 for (
const auto& curQuad : quadLine) {
92 const double integrationElement = intersection.geometry().integrationElement(curQuad.position());
95 const auto uVal = uFunction.evaluate(quadPos);
98 auto neumannValue = neumannBoundaryLoad_(intersection.geometry().global(curQuad.position()), lambda);
100 energy -= neumannValue.dot(uVal) * curQuad.weight() * integrationElement;
106 template <
typename ST>
109 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& dx = std::nullopt)
const {
110 if (not neumannBoundary_ and not neumannBoundaryLoad_)
112 using namespace Dune::DerivativeDirections;
113 using namespace Dune;
114 const auto uFunction = underlying().displacementFunction(par, dx);
116 auto& element = underlying().localView().element();
118 for (
auto&& intersection : intersections(neumannBoundary_->gridView(), element)) {
119 if (not neumannBoundary_->contains(intersection))
123 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), underlying().order());
125 for (
const auto& curQuad : quadLine) {
127 const double intElement = intersection.geometry().integrationElement(curQuad.position());
130 for (
size_t i = 0; i < underlying().numberOfNodes(); ++i) {
131 const auto udCi = uFunction.evaluateDerivative(quadPos, wrt(coeff(i)));
134 auto neumannValue = neumannBoundaryLoad_(intersection.geometry().global(curQuad.position()), lambda);
135 force.template segment<worldDim>(
worldDim * i) -= udCi * neumannValue * curQuad.weight() * intElement;
141 template <
typename ST>
144 const std::optional<std::reference_wrapper<
const Eigen::VectorX<ST>>>& = std::nullopt)
const {}
148 neumannBoundaryLoad_;
149 const BoundaryPatch<GridView>* neumannBoundary_;
152 const auto& underlying()
const {
return static_cast<const FE&
>(*this); }
153 auto& underlying() {
return static_cast<FE&
>(*this); }
164template <
typename GV,
typename F>
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