version 3.8
geomechanics/elastic/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//
13#ifndef DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
14#define DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
15
20
21namespace Dumux {
22
28template<class TypeTag>
29class ElasticLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
30{
32
34 using Element = typename GridView::template Codim<0>::Entity;
35
38 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
39 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
40 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
42 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
44 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
45
46 // class assembling the stress tensor
48
49public:
50 using ParentType::ParentType;
51
60 NumEqVector computeStorage(const Problem& problem,
61 const SubControlVolume& scv,
62 const VolumeVariables& volVars) const
63 {
64 return NumEqVector(0.0);
65 }
66
78 NumEqVector computeFlux(const Problem& problem,
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace& scvf,
83 const ElementFluxVariablesCache& elemFluxVarsCache) const
84 {
85 // obtain force on the face from stress type
86 return StressType::force(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
87 }
88
102 NumEqVector computeSource(const Problem& problem,
103 const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const ElementVolumeVariables& elemVolVars,
106 const SubControlVolume &scv) const
107 {
108 NumEqVector source(0.0);
109
110 // add contributions from volume flux sources
111 source += problem.source(element, fvGeometry, elemVolVars, scv);
112
113 // add contribution from possible point sources
114 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
115
116 // maybe add gravitational acceleration
117 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
118 if (gravity)
119 {
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];
123 }
124
125 return source;
126 }
127};
128
129} // end namespace Dumux
130
131#endif
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
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.