12#if HAVE_DUNE_LOCALFEFUNCTIONS
13# include <dune/fufem/boundarypatch.hh>
14# include <dune/geometry/quadraturerules.hh>
15# include <dune/geometry/type.hh>
16# include <dune/localfefunctions/cachedlocalBasis/cachedlocalBasis.hh>
17# include <dune/localfefunctions/expressions/greenLagrangeStrains.hh>
18# include <dune/localfefunctions/impl/standardLocalFunction.hh>
19# include <dune/localfefunctions/manifolds/realTuple.hh>
43 template <
typename Basis_,
typename Material_,
typename FERequirements_ = FERequirements<>,
bool useEigenRef = false>
53 using Element =
typename LocalView::Element;
60 using LocalBasisType =
decltype(std::declval<LocalView>().tree().child(0).finiteElement().localBasis());
74 template <
typename VolumeLoad = utils::LoadDefault,
typename NeumannBoundaryLoad = utils::LoadDefault>
76 VolumeLoad p_volumeLoad = {},
const BoundaryPatch<GridView>* p_neumannBoundary =
nullptr,
77 NeumannBoundaryLoad p_neumannBoundaryLoad = {})
80 auto& first_child = this->
localView().tree().child(0);
81 const auto& fe = first_child.finiteElement();
82 geo_ = std::make_shared<const Geometry>(this->
localView().element().geometry());
84 order = 2 * (fe.localBasis().order());
85 localBasis = Dune::CachedLocalBasis(fe.localBasis());
86 if constexpr (
requires { this->
localView().element().impl().getQuadratureRule(
order); })
87 if (this->
localView().element().impl().isTrimmed())
88 localBasis.bind(this->localView().element().impl().getQuadratureRule(
order), Dune::bindDerivatives(0, 1));
91 Dune::bindDerivatives(0, 1));
94 Dune::bindDerivatives(0, 1));
96 if constexpr (!std::is_same_v<VolumeLoad, utils::LoadDefault>)
volumeLoad = p_volumeLoad;
97 if constexpr (!std::is_same_v<NeumannBoundaryLoad, utils::LoadDefault>)
101 &&
"If you pass a Neumann boundary you should also pass the function for the Neumann load!");
112 template <
typename ScalarType =
double>
114 const std::optional<
const Eigen::VectorX<ScalarType>>& dx = std::nullopt)
const {
116 auto disp = Ikarus::FEHelper::localSolutionBlockVector<FERequirementType>(d, this->
localView(), dx);
129 template <
typename ScalarType =
double>
131 const std::optional<
const Eigen::VectorX<ScalarType>>& dx = std::nullopt)
const {
144 template <
typename ScalarType,
int strainDim,
bool voigt = true>
146 if constexpr (std::is_same_v<ScalarType, double>)
147 return mat.template tangentModuli<strainType, voigt>(strain);
149 decltype(
auto) matAD =
mat.template rebind<ScalarType>();
150 return matAD.template tangentModuli<strainType, voigt>(strain);
162 template <
typename ScalarType,
int strainDim>
164 if constexpr (std::is_same_v<ScalarType, double>)
165 return mat.template storedEnergy<strainType>(strain);
167 decltype(
auto) matAD =
mat.template rebind<ScalarType>();
168 return matAD.template storedEnergy<strainType>(strain);
181 template <
typename ScalarType,
int strainDim,
bool voigt = true>
182 auto getStress(
const Eigen::Vector<ScalarType, strainDim>& strain)
const {
183 if constexpr (std::is_same_v<ScalarType, double>)
184 return mat.template stresses<strainType, voigt>(strain);
186 decltype(
auto) matAD =
mat.template rebind<ScalarType>();
187 return matAD.template stresses<strainType, voigt>(strain);
208 calculateVectorImpl<double>(par, force);
219 using namespace Dune::DerivativeDirections;
220 using namespace Dune;
222 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
223 const double intElement =
geo_->integrationElement(gp.position()) * gp.weight();
224 const auto EVoigt = (eps.evaluate(gpIndex, on(
gridElement))).eval();
228 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(
gridElement));
230 const auto bopJ = eps.evaluateDerivative(gpIndex, wrt(coeff(j)), on(
gridElement));
231 const auto kgIJ = eps.evaluateDerivative(gpIndex, wrt(coeff(i, j)), along(stresses), on(
gridElement));
232 K.template block<myDim, myDim>(i *
myDim, j *
myDim) += (bopI.transpose() * C * bopJ + kgIJ) * intElement;
247 using namespace Dune::DerivativeDirections;
248 using namespace Dune;
251 const auto H = uFunction.evaluateDerivative(local, Dune::wrt(spatialAll), Dune::on(
gridElement));
252 const auto E = (0.5 * (H.transpose() + H + H.transpose() * H)).eval();
253 const auto EVoigt =
toVoigt(E);
254 auto PK2 =
mat.template stresses<StrainTags::greenLagrangian>(EVoigt);
259 DUNE_THROW(Dune::NotImplemented,
"The requested result type is NOT implemented.");
262 std::shared_ptr<const Geometry>
geo_;
263 Dune::CachedLocalBasis<std::remove_cvref_t<LocalBasisType>>
localBasis;
274 template <
typename ScalarType>
276 = std::nullopt)
const -> ScalarType {
277 using namespace Dune::DerivativeDirections;
278 using namespace Dune;
282 ScalarType energy = 0.0;
284 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
285 const auto EVoigt = (eps.evaluate(gpIndex, on(
gridElement))).eval();
287 energy += internalEnergy *
geo_->integrationElement(gp.position()) * gp.weight();
291 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
292 const auto u = uFunction.evaluate(gpIndex);
293 const Eigen::Vector<double, Traits::worlddim> fExt =
volumeLoad(
geo_->global(gp.position()), lambda);
294 energy -= u.dot(fExt) *
geo_->integrationElement(gp.position()) * gp.weight();
301 const auto& element = this->
localView().element();
302 for (
auto&& intersection : intersections(
neumannBoundary->gridView(), element)) {
305 const auto& quadLine = Dune::QuadratureRules<double, Traits::mydim - 1>::rule(intersection.type(),
order);
307 for (
const auto& curQuad : quadLine) {
310 = intersection.geometryInInside().global(curQuad.position());
312 const double intElement = intersection.geometry().integrationElement(curQuad.position());
315 const auto u = uFunction.evaluate(quadPos);
318 const auto neumannValue =
neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
319 energy -= neumannValue.dot(u) * curQuad.weight() * intElement;
325 template <
typename ScalarType>
327 const std::optional<
const Eigen::VectorX<ScalarType>>& dx = std::nullopt)
const {
328 using namespace Dune::DerivativeDirections;
329 using namespace Dune;
335 for (
const auto& [gpIndex, gp] : eps.viewOverIntegrationPoints()) {
336 const double intElement =
geo_->integrationElement(gp.position()) * gp.weight();
337 const auto EVoigt = (eps.evaluate(gpIndex, on(
gridElement))).eval();
340 const auto bopI = eps.evaluateDerivative(gpIndex, wrt(coeff(i)), on(
gridElement));
341 force.template segment<myDim>(
myDim * i) += bopI.transpose() * stresses * intElement;
347 for (
const auto& [gpIndex, gp] : uFunction.viewOverIntegrationPoints()) {
348 const double intElement =
geo_->integrationElement(gp.position()) * gp.weight();
349 const Eigen::Vector<double, Traits::worlddim> fExt =
volumeLoad(
geo_->global(gp.position()), lambda);
351 const auto udCi = uFunction.evaluateDerivative(gpIndex, wrt(coeff(i)));
352 force.template segment<myDim>(
myDim * i) -= udCi * fExt * intElement;
359 const auto& element = this->
localView().element();
360 for (
auto&& intersection : intersections(
neumannBoundary->gridView(), element)) {
364 const auto& quadLine = Dune::QuadratureRules<double, myDim - 1>::rule(intersection.type(),
order);
366 for (
const auto& curQuad : quadLine) {
369 const double intElement = intersection.geometry().integrationElement(curQuad.position());
373 const auto udCi = uFunction.evaluateDerivative(quadPos, wrt(coeff(i)));
376 const auto neumannValue =
neumannBoundaryLoad(intersection.geometry().global(curQuad.position()), lambda);
377 force.template segment<myDim>(
myDim * i) -= udCi * neumannValue * curQuad.weight() * intElement;
386# error NonLinearElastic depends on dune-localfefunctions, which is not included
Collection of fallback default functions.
Helper for transform between Dune linear algebra types and Eigen.
Helper for the autodiff library.
Contains the PowerBasisFE class, which works with a power basis in FlatInterLeaved elements.
Material property functions and conversion utilities.
FETraits template structure for finite element traits.
Definition of the LinearElastic class for finite element mechanics computations.
Definition of several material related enums.
constexpr Eigen::Index toVoigt(Eigen::Index i, Eigen::Index j) noexcept
Converts 2D indices to Voigt notation index.
Definition: tensorutils.hh:167
Definition: simpleassemblers.hh:21
Definition: resultevaluators.hh:17
PowerBasisFE class for working with a power basis in FlatInterLeaved elements.
Definition: powerbasisfe.hh:23
const GridElement & gridElement() const
Get the grid element associated with the local view.
Definition: powerbasisfe.hh:78
const LocalView & localView() const
Get the const reference to the local view.
Definition: powerbasisfe.hh:84
Class representing a map of result types to result arrays.
Definition: ferequirements.hh:342
void insertOrAssignResult(ResultType &&resultType, const ResultArray &resultArray)
Insert or assign a result to the map.
Definition: ferequirements.hh:354
Class representing the requirements for obtaining specific results.
Definition: ferequirements.hh:405
const FERequirements & getFERequirements() const
Get the associated finite element requirements.
Definition: ferequirements.hh:535
bool isResultRequested(ResultType &&key) const
Check if a specific result type is requested.
Definition: ferequirements.hh:446
NonLinearElastic class represents a non-linear elastic finite element.
Definition: finiteelements/mechanics/nonlinearelastic.hh:44
static constexpr int myDim
Definition: finiteelements/mechanics/nonlinearelastic.hh:57
typename LocalView::Element Element
Definition: finiteelements/mechanics/nonlinearelastic.hh:53
static constexpr auto strainType
Definition: finiteelements/mechanics/nonlinearelastic.hh:59
std::shared_ptr< const Geometry > geo_
Definition: finiteelements/mechanics/nonlinearelastic.hh:262
Basis_ Basis
Definition: finiteelements/mechanics/nonlinearelastic.hh:46
PowerBasisFE< FlatBasis > BasePowerFE
Definition: finiteelements/mechanics/nonlinearelastic.hh:49
void calculateMatrix(const FERequirementType &par, typename Traits::template MatrixType<> K) const
Calculate the matrix associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:218
auto getInternalEnergy(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the internal energy for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:163
double calculateScalar(const FERequirementType &par) const
Calculate the scalar value associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:198
typename FlatBasis::GridView GridView
Definition: finiteelements/mechanics/nonlinearelastic.hh:55
std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> volumeLoad
Definition: finiteelements/mechanics/nonlinearelastic.hh:265
auto displacementFunction(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Get the displacement function for the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:113
Dune::CachedLocalBasis< std::remove_cvref_t< LocalBasisType > > localBasis
Definition: finiteelements/mechanics/nonlinearelastic.hh:263
typename Element::Geometry Geometry
Definition: finiteelements/mechanics/nonlinearelastic.hh:54
auto calculateScalarImpl(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const -> ScalarType
Definition: finiteelements/mechanics/nonlinearelastic.hh:275
FERequirements_ FERequirementType
Definition: finiteelements/mechanics/nonlinearelastic.hh:50
void calculateVectorImpl(const FERequirementType &par, typename Traits::template VectorType< ScalarType > force, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
Definition: finiteelements/mechanics/nonlinearelastic.hh:326
Material mat
Definition: finiteelements/mechanics/nonlinearelastic.hh:269
size_t numberOfNodes
Definition: finiteelements/mechanics/nonlinearelastic.hh:270
Material_ Material
Definition: finiteelements/mechanics/nonlinearelastic.hh:47
decltype(std::declval< LocalView >().tree().child(0).finiteElement().localBasis()) LocalBasisType
Definition: finiteelements/mechanics/nonlinearelastic.hh:60
std::function< Eigen::Vector< double, worldDim >(const Dune::FieldVector< double, worldDim > &, const double &)> neumannBoundaryLoad
Definition: finiteelements/mechanics/nonlinearelastic.hh:267
const BoundaryPatch< GridView > * neumannBoundary
Definition: finiteelements/mechanics/nonlinearelastic.hh:268
auto getStress(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the stress for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:182
typename Basis::FlatBasis FlatBasis
Definition: finiteelements/mechanics/nonlinearelastic.hh:48
static constexpr int worldDim
Definition: finiteelements/mechanics/nonlinearelastic.hh:58
auto materialTangent(const Eigen::Vector< ScalarType, strainDim > &strain) const
Get the material tangent for the given strain.
Definition: finiteelements/mechanics/nonlinearelastic.hh:145
auto strainFunction(const FERequirementType &par, const std::optional< const Eigen::VectorX< ScalarType > > &dx=std::nullopt) const
The strain function for the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:130
typename FlatBasis::LocalView LocalView
Definition: finiteelements/mechanics/nonlinearelastic.hh:52
int order
Definition: finiteelements/mechanics/nonlinearelastic.hh:271
void calculateAt(const ResultRequirementsType &req, const Dune::FieldVector< double, Traits::mydim > &local, ResultTypeMap< double > &result) const
Calculate specified results at a given local position.
Definition: finiteelements/mechanics/nonlinearelastic.hh:245
NonLinearElastic(const Basis &globalBasis, const typename LocalView::Element &element, const Material &p_mat, VolumeLoad p_volumeLoad={}, const BoundaryPatch< GridView > *p_neumannBoundary=nullptr, NeumannBoundaryLoad p_neumannBoundaryLoad={})
Constructor for the NonLinearElastic class.
Definition: finiteelements/mechanics/nonlinearelastic.hh:75
void calculateVector(const FERequirementType &par, typename Traits::template VectorType<> force) const
Calculate the vector associated with the given FERequirementType.
Definition: finiteelements/mechanics/nonlinearelastic.hh:207
Traits for handling local views.see https://en.wikipedia.org/wiki/Lam%C3%A9_parameters.
Definition: physicshelper.hh:65
static constexpr int worlddim
Dimension of the world space.
Definition: physicshelper.hh:68
static constexpr int mydim
Dimension of the geometry.
Definition: physicshelper.hh:71
Definition: resultevaluators.hh:20
decltype(Dune::Functions::DefaultGlobalBasis(Ikarus::flatPreBasis(std::declval< PreBasis >()))) FlatBasis
The type of the flattened basis.
Definition: utils/basis.hh:36