6#include <dune/fufem/boundarypatch.hh>
21 template <
typename DisplacementBasedElement,
typename Traits>
27 static constexpr int myDim = Traits::mydim;
28 static constexpr int worldDim = Traits::worlddim;
37 template <
typename NeumannBoundaryLoad>
38 explicit Traction(
const BoundaryPatch<GridView>* p_neumannBoundary, NeumannBoundaryLoad p_neumannBoundaryLoad)
39 : neumannBoundary{p_neumannBoundary} {
40 if constexpr (!std::is_same_v<NeumannBoundaryLoad, utils::LoadDefault>)
41 neumannBoundaryLoad = p_neumannBoundaryLoad;
43 assert(((not p_neumannBoundary and not neumannBoundaryLoad) or (p_neumannBoundary and neumannBoundaryLoad))
44 &&
"If you pass a Neumann boundary you should also pass the function for the Neumann load!");
65 calculateVectorImpl<double>(req, force);
75 calculateMatrixImpl<double>(req, K);
79 template <
typename ScalarType>
81 = std::nullopt)
const -> ScalarType {
82 if (not neumannBoundary and not neumannBoundaryLoad)
return 0.0;
83 ScalarType energy = 0.0;
84 const auto uFunction = dbElement().displacementFunction(par, dx);
86 auto& element = dbElement().localView().element();
88 for (
auto&& intersection : intersections(neumannBoundary->gridView(), element)) {
89 if (not neumannBoundary->contains(intersection))
continue;
91 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), dbElement().order());
93 for (
const auto& curQuad : quadLine) {
97 const double integrationElement = intersection.geometry().integrationElement(curQuad.position());
100 const auto uVal = uFunction.evaluate(quadPos);
103 auto neumannValue = neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
105 energy -= neumannValue.dot(uVal) * curQuad.weight() * integrationElement;
111 template <
typename ScalarType>
113 const std::optional<
const Eigen::VectorX<ScalarType>> dx = std::nullopt)
const {
114 if (not neumannBoundary and not neumannBoundaryLoad)
return;
115 using namespace Dune::DerivativeDirections;
116 using namespace Dune;
117 const auto uFunction = dbElement().displacementFunction(par, dx);
119 auto& element = dbElement().localView().element();
121 for (
auto&& intersection : intersections(neumannBoundary->gridView(), element)) {
122 if (not neumannBoundary->contains(intersection))
continue;
125 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), dbElement().order());
127 for (
const auto& curQuad : quadLine) {
129 const double intElement = intersection.geometry().integrationElement(curQuad.position());
132 for (
size_t i = 0; i < dbElement().numberOfNodes(); ++i) {
133 const auto udCi = uFunction.evaluateDerivative(quadPos, wrt(coeff(i)));
136 auto neumannValue = neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
137 force.template segment<worldDim>(
worldDim * i) -= udCi * neumannValue * curQuad.weight() * intElement;
143 template <
typename ScalarType>
145 const std::optional<
const Eigen::VectorX<ScalarType>>& dx = std::nullopt)
const {}
150 const BoundaryPatch<GridView>* neumannBoundary;
157 constexpr const DisplacementBasedElement& dbElement()
const {
158 return static_cast<DisplacementBasedElement const&
>(*this);
Collection of fallback default functions.
Definition of the LinearElastic class for finite element mechanics computations.
Definition: simpleassemblers.hh:21
Definition: resultevaluators.hh:17
Traction class represents distributed traction load that can be applied.
Definition: traction.hh:22
typename Traits::GridView GridView
Definition: traction.hh:26
void calculateVectorImpl(const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > dx=std::nullopt) const
Definition: traction.hh:112
static constexpr int worldDim
Definition: traction.hh:28
static constexpr int myDim
Definition: traction.hh:27
auto calculateScalarImpl(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType
Definition: traction.hh:80
void calculateVector(const FERequirementType &req, typename Traits::template VectorType<> force) const
Calculate the vector associated with the given FERequirementType.
Definition: traction.hh:64
typename Traits::LocalView LocalView
Definition: traction.hh:25
void calculateMatrixImpl(const FERequirementType &par, typename Traits::template MatrixType<> K, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Definition: traction.hh:144
Traction(const BoundaryPatch< GridView > *p_neumannBoundary, NeumannBoundaryLoad p_neumannBoundaryLoad)
Constructor for the Loads class.
Definition: traction.hh:38
typename Traits::FERequirementType FERequirementType
Definition: traction.hh:24
double calculateScalar(const FERequirementType &req) const
Calculate the scalar value.
Definition: traction.hh:55
void calculateMatrix(const FERequirementType &req, typename Traits::template MatrixType<> K) const
Calculate the matrix associated with the given FERequirementType.
Definition: traction.hh:74
Definition: resultevaluators.hh:20