3.6-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
32
33namespace Dumux {
34
40template<class TypeTag>
41class ElasticLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
42{
44
46 using Element = typename GridView::template Codim<0>::Entity;
47
50 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
51 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
57
58 // class assembling the stress tensor
60
61public:
62 using ParentType::ParentType;
63
72 NumEqVector computeStorage(const Problem& problem,
73 const SubControlVolume& scv,
74 const VolumeVariables& volVars) const
75 {
76 return NumEqVector(0.0);
77 }
78
90 NumEqVector computeFlux(const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const SubControlVolumeFace& scvf,
95 const ElementFluxVariablesCache& elemFluxVarsCache) const
96 {
97 // obtain force on the face from stress type
98 return StressType::force(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
99 }
100
114 NumEqVector computeSource(const Problem& problem,
115 const Element& element,
116 const FVElementGeometry& fvGeometry,
117 const ElementVolumeVariables& elemVolVars,
118 const SubControlVolume &scv) const
119 {
120 NumEqVector source(0.0);
121
122 // add contributions from volume flux sources
123 source += problem.source(element, fvGeometry, elemVolVars, scv);
124
125 // add contribution from possible point sources
126 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
127
128 // maybe add gravitational acceleration
129 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
130 if (gravity)
131 {
132 const auto& g = problem.spatialParams().gravity(scv.center());
133 for (int dir = 0; dir < GridView::dimensionworld; ++dir)
134 source[Indices::momentum(dir)] += elemVolVars[scv].solidDensity()*g[dir];
135 }
136
137 return source;
138 }
139};
140
141} // end namespace Dumux
142
143#endif
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Element-wise calculation of the local residual for problems using the elastic model considering linea...
Definition: geomechanics/elastic/localresidual.hh:42
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:90
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:72
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:114
Declares all properties used in Dumux.