3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
26#define DUMUX_GEOMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH
27
31
32namespace Dumux {
33
39template<class TypeTag>
40class ElasticLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
41{
43
45 using Element = typename GridView::template Codim<0>::Entity;
46
49 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
50 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
53 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
56
57 // class assembling the stress tensor
59
60public:
61 using ParentType::ParentType;
62
71 NumEqVector computeStorage(const Problem& problem,
72 const SubControlVolume& scv,
73 const VolumeVariables& volVars) const
74 {
75 return NumEqVector(0.0);
76 }
77
89 NumEqVector computeFlux(const Problem& problem,
90 const Element& element,
91 const FVElementGeometry& fvGeometry,
92 const ElementVolumeVariables& elemVolVars,
93 const SubControlVolumeFace& scvf,
94 const ElementFluxVariablesCache& elemFluxVarsCache) const
95 {
96 // obtain force on the face from stress type
97 return StressType::force(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
98 }
99
113 NumEqVector computeSource(const Problem& problem,
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& elemVolVars,
117 const SubControlVolume &scv) const
118 {
119 NumEqVector source(0.0);
120
121 // add contributions from volume flux sources
122 source += problem.source(element, fvGeometry, elemVolVars, scv);
123
124 // add contribution from possible point sources
125 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
126
127 // maybe add gravitational acceleration
128 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
129 if (gravity)
130 {
131 const auto& g = problem.spatialParams().gravity(scv.center());
132 for (int dir = 0; dir < GridView::dimensionworld; ++dir)
133 source[Indices::momentum(dir)] += elemVolVars[scv].solidDensity()*g[dir];
134 }
135
136 return source;
137 }
138};
139
140} // end namespace Dumux
141
142#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Element-wise calculation of the local residual for problems using the elastic model considering linea...
Definition: geomechanics/elastic/localresidual.hh:41
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:89
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:71
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:113
Declares all properties used in Dumux.