3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/box/darcyslaw.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 *****************************************************************************/
28#ifndef DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
30
31#include <dumux/common/math.hh>
36
37namespace Dumux {
38
39// forward declaration
40template<class TypeTag, DiscretizationMethod discMethod>
41class DarcysLawImplementation;
42
43// forward declaration
44template<class Scalar, class GridGeometry>
45class BoxDarcysLaw;
46
51template<class TypeTag>
53: public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
54{ };
55
62template<class Scalar, class GridGeometry>
64{
65 using FVElementGeometry = typename GridGeometry::LocalView;
66 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
67 using Extrusion = Extrusion_t<GridGeometry>;
68 using GridView = typename GridGeometry::GridView;
69 using Element = typename GridView::template Codim<0>::Entity;
70
71 enum { dim = GridView::dimension};
72 enum { dimWorld = GridView::dimensionworld};
73
74public:
75
86 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
87 static Scalar flux(const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
92 const int phaseIdx,
93 const ElementFluxVarsCache& elemFluxVarCache)
94 {
95 const auto& fluxVarCache = elemFluxVarCache[scvf];
96 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
97 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
98 const auto& insideVolVars = elemVolVars[insideScv];
99 const auto& outsideVolVars = elemVolVars[outsideScv];
100
101 auto insideK = insideVolVars.permeability();
102 auto outsideK = outsideVolVars.permeability();
103
104 // scale with correct extrusion factor
105 insideK *= insideVolVars.extrusionFactor();
106 outsideK *= outsideVolVars.extrusionFactor();
107
108 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
109 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
110
111 const auto& shapeValues = fluxVarCache.shapeValues();
112
113 // evaluate gradP - rho*g at integration point
114 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
115 Scalar rho(0.0);
116 for (auto&& scv : scvs(fvGeometry))
117 {
118 const auto& volVars = elemVolVars[scv];
119
120 if (enableGravity)
121 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
122
123 // the global shape function gradient
124 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
125 }
126
127 if (enableGravity)
128 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
129
130 // apply the permeability and return the flux
131 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(scvf);
132 }
133
134 // compute transmissibilities ti for analytical jacobians
135 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
136 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const FluxVarCache& fluxVarCache)
142 {
143 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
144 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
145 const auto& insideVolVars = elemVolVars[insideScv];
146 const auto& outsideVolVars = elemVolVars[outsideScv];
147
148 auto insideK = insideVolVars.permeability();
149 auto outsideK = outsideVolVars.permeability();
150
151 // scale with correct extrusion factor
152 insideK *= insideVolVars.extrusionFactor();
153 outsideK *= outsideVolVars.extrusionFactor();
154
155 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
156
157 std::vector<Scalar> ti(fvGeometry.numScv());
158 for (const auto& scv : scvs(fvGeometry))
159 ti[scv.indexInElement()] =
160 -1.0*Extrusion::area(scvf)*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
161
162 return ti;
163 }
164};
165
166} // end namespace Dumux
167
168#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:64
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache)
Definition: flux/box/darcyslaw.hh:136
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Returns the advective flux of a fluid phase across the given sub-control volume face.
Definition: flux/box/darcyslaw.hh:87
Declares all properties used in Dumux.