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>
29template<
class TypeTag>
31:
public GetPropType<TypeTag, Properties::BaseLocalResidual>
41 using FVElementGeometry =
typename GridGeometry::LocalView;
42 using SubControlVolume =
typename GridGeometry::SubControlVolume;
43 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
44 using GridView =
typename GridGeometry::GridView;
45 using Element =
typename GridView::template Codim<0>::Entity;
48 using Indices =
typename ModelTraits::Indices;
49 static constexpr int dimWorld = GridView::dimensionworld;
51 using ParentType::ParentType;
58 const SubControlVolume& scv,
59 const VolumeVariables& volVars)
const
62 return NumEqVector(0.0);
70 const Element& element,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
73 const SubControlVolumeFace& scvf,
74 const ElementFluxVariablesCache& elemFluxVarsCache)
const
79 const auto& fluxVarCache = elemFluxVarsCache[scvf];
80 Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0);
81 for (
const auto& scv : scvs(fvGeometry))
83 const auto& volVars = elemVolVars[scv];
84 for (
int dir = 0; dir < dimWorld; ++dir)
85 F[dir].axpy(volVars.displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
88 for (
int dir = 0; dir < dimWorld; ++dir)
95 NumEqVector flux(0.0);
98 const auto P = problem.firstPiolaKirchhoffStressTensor(F);
99 P.mv(scvf.unitOuterNormal(), flux);
100 flux *= -scvf.area();
105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const SubControlVolume& scv)
const
111 return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
Local residual for the hyperelastic model.
Definition: geomechanics/hyperelastic/localresidual.hh:32
typename ParentType::ElementResidualVector ElementResidualVector
Definition: geomechanics/hyperelastic/localresidual.hh:52
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate the rate of change of all conserved quantities.
Definition: geomechanics/hyperelastic/localresidual.hh:57
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Definition: geomechanics/hyperelastic/localresidual.hh:104
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:69
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
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.