13#ifndef DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
14#define DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
28template<
class TypeTag>
34 using Element =
typename GridView::template Codim<0>::Entity;
39 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
40 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
44 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
50 using ParentType::ParentType;
61 const SubControlVolume& scv,
62 const VolumeVariables& volVars)
const
64 return NumEqVector(0.0);
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace& scvf,
83 const ElementFluxVariablesCache& elemFluxVarsCache)
const
86 return StressType::force(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
103 const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const ElementVolumeVariables& elemVolVars,
106 const SubControlVolume &scv)
const
108 NumEqVector source(0.0);
111 source += problem.source(element, fvGeometry, elemVolVars, scv);
114 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
117 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
120 const auto& g = problem.spatialParams().gravity(scv.center());
121 for (
int dir = 0; dir < GridView::dimensionworld; ++dir)
122 source[Indices::momentum(dir)] += elemVolVars[scv].solidDensity()*g[dir];
Element-wise calculation of the local residual for problems using the elastic model considering linea...
Definition: geomechanics/elastic/localresidual.hh:30
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the force in all grid directions acting on a face of a sub-control volume.
Definition: geomechanics/elastic/localresidual.hh:78
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
For the elastic model the storage term is zero since we neglect inertia forces.
Definition: geomechanics/elastic/localresidual.hh:60
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: geomechanics/elastic/localresidual.hh:102
Defines all properties used in Dumux.
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
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.