version 0.4
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
10
11namespace Ikarus {
12
21 template <typename DisplacementBasedElement, typename Traits>
22 class Traction {
23 public:
24 using FERequirementType = typename Traits::FERequirementType;
25 using LocalView = typename Traits::LocalView;
26 using GridView = typename Traits::GridView;
27 static constexpr int myDim = Traits::mydim;
28 static constexpr int worldDim = Traits::worlddim;
29
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;
42
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!");
45 }
46
55 double calculateScalar(const FERequirementType& req) const { return calculateScalarImpl<double>(req); }
56
64 void calculateVector(const FERequirementType& req, typename Traits::template VectorType<> force) const {
65 calculateVectorImpl<double>(req, force);
66 }
74 void calculateMatrix(const FERequirementType& req, typename Traits::template MatrixType<> K) const {
75 calculateMatrixImpl<double>(req, K);
76 }
77
78 protected:
79 template <typename ScalarType>
80 auto calculateScalarImpl(const FERequirementType& par, const std::optional<const Eigen::VectorX<ScalarType>>& dx
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);
85 const auto& lambda = par.getParameter(Ikarus::FEParameter::loadfactor);
86 auto& element = dbElement().localView().element();
87
88 for (auto&& intersection : intersections(neumannBoundary->gridView(), element)) {
89 if (not neumannBoundary->contains(intersection)) continue;
90
91 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), dbElement().order());
92
93 for (const auto& curQuad : quadLine) {
94 // Local position of the quadrature point
95 const Dune::FieldVector<double, myDim>& quadPos = intersection.geometryInInside().global(curQuad.position());
96
97 const double integrationElement = intersection.geometry().integrationElement(curQuad.position());
98
99 // The value of the local function
100 const auto uVal = uFunction.evaluate(quadPos);
101
102 // Value of the Neumann data at the current position
103 auto neumannValue = neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
104
105 energy -= neumannValue.dot(uVal) * curQuad.weight() * integrationElement;
106 }
107 }
108 return energy;
109 }
110
111 template <typename ScalarType>
112 void calculateVectorImpl(const FERequirementType& par, typename Traits::template VectorType<ScalarType> force,
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);
118 const auto& lambda = par.getParameter(Ikarus::FEParameter::loadfactor);
119 auto& element = dbElement().localView().element();
120
121 for (auto&& intersection : intersections(neumannBoundary->gridView(), element)) {
122 if (not neumannBoundary->contains(intersection)) continue;
123
125 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(), dbElement().order());
126
127 for (const auto& curQuad : quadLine) {
128 const Dune::FieldVector<double, myDim>& quadPos = intersection.geometryInInside().global(curQuad.position());
129 const double intElement = intersection.geometry().integrationElement(curQuad.position());
130
132 for (size_t i = 0; i < dbElement().numberOfNodes(); ++i) {
133 const auto udCi = uFunction.evaluateDerivative(quadPos, wrt(coeff(i)));
134
136 auto neumannValue = neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
137 force.template segment<worldDim>(worldDim * i) -= udCi * neumannValue * curQuad.weight() * intElement;
138 }
139 }
140 }
141 }
142
143 template <typename ScalarType>
144 void calculateMatrixImpl(const FERequirementType& par, typename Traits::template MatrixType<> K,
145 const std::optional<const Eigen::VectorX<ScalarType>>& dx = std::nullopt) const {}
146
147 private:
148 std::function<Eigen::Vector<double, worldDim>(const Dune::FieldVector<double, worldDim>&, const double&)>
149 neumannBoundaryLoad;
150 const BoundaryPatch<GridView>* neumannBoundary;
151
157 constexpr const DisplacementBasedElement& dbElement() const {
158 return static_cast<DisplacementBasedElement const&>(*this);
159 }
160 };
161} // namespace Ikarus
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