version 3.10-dev
geomechanics/hyperelastic/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH
13#define DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH
14
15#include <dune/common/fvector.hh>
16#include <dune/common/fmatrix.hh>
17#include <dune/common/std/type_traits.hh>
18#include <dune/common/exceptions.hh>
19
20#include <dumux/common/math.hh>
24
25namespace Dumux {
26
27#ifndef DOXYGEN
28namespace Detail::Hyperelastic {
29// helper struct detecting if the user-defined problem class has a solidDensity method
30template <typename T, typename ...Ts>
31using SolidDensityDetector = decltype(std::declval<T>().solidDensity(std::declval<Ts>()...));
32
33template<class T, typename ...Args>
34static constexpr bool hasSolidDensity()
35{ return Dune::Std::is_detected<SolidDensityDetector, T, Args...>::value; }
36
37// helper struct detecting if the user-defined problem class has a acceleration method
38template <typename T, typename ...Ts>
39using AccelerationDetector = decltype(std::declval<T>().acceleration(std::declval<Ts>()...));
40
41template<class T, typename ...Args>
42static constexpr bool hasAcceleration()
43{ return Dune::Std::is_detected<AccelerationDetector, T, Args...>::value; }
44} // end namespace Detail::Hyperelastic
45#endif
46
51template<class TypeTag>
53: public GetPropType<TypeTag, Properties::BaseLocalResidual>
54{
60 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
63 using FVElementGeometry = typename GridGeometry::LocalView;
64 using SubControlVolume = typename GridGeometry::SubControlVolume;
65 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
66 using GridView = typename GridGeometry::GridView;
67 using Element = typename GridView::template Codim<0>::Entity;
68 using Extrusion = Extrusion_t<GridGeometry>;
70 using Indices = typename ModelTraits::Indices;
71 static constexpr int dimWorld = GridView::dimensionworld;
72public:
73 using ParentType::ParentType;
74 using ElementResidualVector = typename ParentType::ElementResidualVector;
75
76 using ParentType::evalStorage;
77
82 const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& prevElemVolVars,
86 const ElementVolumeVariables& curElemVolVars,
87 const SubControlVolume& scv) const
88 {
89 const auto& curVolVars = curElemVolVars[scv];
90
91 if constexpr (Detail::Hyperelastic::hasSolidDensity<Problem, const Element&, const SubControlVolume&>()
92 && Detail::Hyperelastic::hasAcceleration<Problem, const Element&, const SubControlVolume&, Scalar, decltype(curVolVars.displacement())>())
93 {
94 const auto& curVolVars = curElemVolVars[scv];
95 NumEqVector storage = problem.solidDensity(element, scv)*problem.acceleration(
96 element, scv, this->timeLoop().timeStepSize(), curVolVars.displacement()
97 );
98 storage *= curVolVars.extrusionFactor()*Extrusion::volume(fvGeometry, scv);
99 residual[scv.localDofIndex()] += storage;
100 }
101 else
102 DUNE_THROW(Dune::InvalidStateException,
103 "Problem class must implement solidDensity and acceleration"
104 " methods to evaluate storage term"
105 );
106 }
107
112 NumEqVector computeFlux(const Problem& problem,
113 const Element& element,
114 const FVElementGeometry& fvGeometry,
115 const ElementVolumeVariables& elemVolVars,
116 const SubControlVolumeFace& scvf,
117 const ElementFluxVariablesCache& elemFluxVarsCache) const
118 {
119 // the deformation gradient F = grad(d) + I
120 const auto F = [&]()
121 {
122 const auto& fluxVarCache = elemFluxVarsCache[scvf];
123 Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0);
124 for (const auto& scv : scvs(fvGeometry))
125 {
126 const auto& volVars = elemVolVars[scv];
127 for (int dir = 0; dir < dimWorld; ++dir)
128 F[dir].axpy(volVars.displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
129 }
130
131 for (int dir = 0; dir < dimWorld; ++dir)
132 F[dir][dir] += 1;
133
134 return F;
135 }();
136
137 // numerical flux ==> -P*n*dA
138 NumEqVector flux(0.0);
139
140 // problem implements the constitutive law P=P(F)
141 const auto P = problem.firstPiolaKirchhoffStressTensor(F);
142 P.mv(scvf.unitOuterNormal(), flux);
143 flux *= -scvf.area();
144 return flux;
145 }
146
147 NumEqVector computeSource(const Problem& problem,
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const ElementVolumeVariables& elemVolVars,
151 const SubControlVolume& scv) const
152 {
153 // loading and body forces are implemented in the problem
154 return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
155 }
156};
157
158} // end namespace Dumux
159
160#endif
Local residual for the hyperelastic model.
Definition: geomechanics/hyperelastic/localresidual.hh:54
typename ParentType::ElementResidualVector ElementResidualVector
Definition: geomechanics/hyperelastic/localresidual.hh:74
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Definition: geomechanics/hyperelastic/localresidual.hh:147
void evalStorage(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Evaluate the storage contribution to the residual: rho*d^2u/dt^2.
Definition: geomechanics/hyperelastic/localresidual.hh:81
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the fluxes over a face of a sub control volume flux: -div(P) where P is the 1st Piola-Kircho...
Definition: geomechanics/hyperelastic/localresidual.hh:112
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.