12#ifndef DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH
13#define DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH
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>
28namespace Detail::Hyperelastic {
30template <
typename T,
typename ...Ts>
31using SolidDensityDetector =
decltype(std::declval<T>().solidDensity(std::declval<Ts>()...));
33template<
class T,
typename ...Args>
34static constexpr bool hasSolidDensity()
35{
return Dune::Std::is_detected<SolidDensityDetector, T, Args...>::value; }
38template <
typename T,
typename ...Ts>
39using AccelerationDetector =
decltype(std::declval<T>().acceleration(std::declval<Ts>()...));
41template<
class T,
typename ...Args>
42static constexpr bool hasAcceleration()
43{
return Dune::Std::is_detected<AccelerationDetector, T, Args...>::value; }
51template<
class TypeTag>
53:
public GetPropType<TypeTag, Properties::BaseLocalResidual>
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;
70 using Indices =
typename ModelTraits::Indices;
71 static constexpr int dimWorld = GridView::dimensionworld;
73 using ParentType::ParentType;
76 using ParentType::evalStorage;
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
89 const auto& curVolVars = curElemVolVars[scv];
91 if constexpr (Detail::Hyperelastic::hasSolidDensity<Problem, const Element&, const SubControlVolume&>()
92 && Detail::Hyperelastic::hasAcceleration<Problem,
const Element&,
const SubControlVolume&, Scalar,
decltype(curVolVars.displacement())>())
94 const auto& curVolVars = curElemVolVars[scv];
95 NumEqVector storage = problem.solidDensity(element, scv)*problem.acceleration(
96 element, scv, this->timeLoop().timeStepSize(), curVolVars.displacement()
99 residual[scv.localDofIndex()] += storage;
102 DUNE_THROW(Dune::InvalidStateException,
103 "Problem class must implement solidDensity and acceleration"
104 " methods to evaluate storage term"
113 const Element& element,
114 const FVElementGeometry& fvGeometry,
115 const ElementVolumeVariables& elemVolVars,
116 const SubControlVolumeFace& scvf,
117 const ElementFluxVariablesCache& elemFluxVarsCache)
const
122 const auto& fluxVarCache = elemFluxVarsCache[scvf];
123 Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0);
124 for (
const auto& scv : scvs(fvGeometry))
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()));
131 for (
int dir = 0; dir < dimWorld; ++dir)
138 NumEqVector flux(0.0);
141 const auto P = problem.firstPiolaKirchhoffStressTensor(F);
142 P.mv(scvf.unitOuterNormal(), flux);
143 flux *= -scvf.area();
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const ElementVolumeVariables& elemVolVars,
151 const SubControlVolume& scv)
const
154 return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
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.
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.