version 3.8
geomechanics/poroelastic/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_PRORELASTIC_LOCAL_RESIDUAL_HH
14#define DUMUX_PRORELASTIC_LOCAL_RESIDUAL_HH
15
19
20namespace Dumux {
21
27template<class TypeTag>
29{
30 using ParentType = ElasticLocalResidual<TypeTag>;
31
33 using Element = typename GridView::template Codim<0>::Entity;
34
37 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
38 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
39 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
41 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
42 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
43 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
44
45public:
46 using ParentType::ParentType;
47
61 NumEqVector computeSource(const Problem& problem,
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolume &scv) const
66 {
67 NumEqVector source(0.0);
68
69 // add contributions from volume flux sources
70 source += problem.source(element, fvGeometry, elemVolVars, scv);
71
72 // add contribution from possible point sources
73 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
74
75 // maybe add gravitational acceleration
76 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
77 if (gravity)
78 {
79 // compute average density
80 const auto& vv = elemVolVars[scv];
81 const auto phi = vv.porosity();
82 const auto rhoFluid = problem.spatialParams().effectiveFluidDensity(element, scv);
83 const auto rhoAverage = phi*rhoFluid + (1.0 - phi)*vv.solidDensity();
84
85 // add body force
86 const auto& g = problem.spatialParams().gravity(scv.center());
87 for (int dir = 0; dir < GridView::dimensionworld; ++dir)
88 source[ Indices::momentum(dir) ] += rhoAverage*g[dir];
89 }
90
91 return source;
92 }
93};
94
95} // end namespace Dumux
96
97#endif
Element-wise calculation of the local residual for problems using the elastic model considering linea...
Definition: geomechanics/elastic/localresidual.hh:30
Element-wise calculation of the local residual for problems using the poroelastic model.
Definition: geomechanics/poroelastic/localresidual.hh:29
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/poroelastic/localresidual.hh:61
Defines all properties used in Dumux.
Element-wise calculation of the local residual for problems using the elastic model considering linea...
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.