3.1-git
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>
35
36namespace Dumux {
37
38// forward declaration
39template<class TypeTag, DiscretizationMethod discMethod>
40class DarcysLawImplementation;
41
42// forward declaration
43template<class Scalar, class GridGeometry>
44class BoxDarcysLaw;
45
50template<class TypeTag>
52: public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
53{ };
54
61template<class Scalar, class GridGeometry>
63{
64 using FVElementGeometry = typename GridGeometry::LocalView;
65 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
66 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
67 using GridView = typename GridGeometry::GridView;
68 using Element = typename GridView::template Codim<0>::Entity;
69 using CoordScalar = typename GridView::ctype;
70
71 enum { dim = GridView::dimension};
72 enum { dimWorld = GridView::dimensionworld};
73
74 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
75
76public:
77
78 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
79 static Scalar flux(const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const SubControlVolumeFace& scvf,
84 const int phaseIdx,
85 const ElementFluxVarsCache& elemFluxVarCache)
86 {
87 const auto& fluxVarCache = elemFluxVarCache[scvf];
88 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
89 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
90 const auto& insideVolVars = elemVolVars[insideScv];
91 const auto& outsideVolVars = elemVolVars[outsideScv];
92
93 auto insideK = insideVolVars.permeability();
94 auto outsideK = outsideVolVars.permeability();
95
96 // scale with correct extrusion factor
97 insideK *= insideVolVars.extrusionFactor();
98 outsideK *= outsideVolVars.extrusionFactor();
99
100 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
101 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
102
103 const auto& shapeValues = fluxVarCache.shapeValues();
104
105 // evaluate gradP - rho*g at integration point
106 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
107 Scalar rho(0.0);
108 for (auto&& scv : scvs(fvGeometry))
109 {
110 const auto& volVars = elemVolVars[scv];
111
112 if (enableGravity)
113 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
114
115 // the global shape function gradient
116 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
117 }
118
119 if (enableGravity)
120 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
121
122 // apply the permeability and return the flux
123 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*scvf.area();
124 }
125
126 // compute transmissibilities ti for analytical jacobians
127 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
128 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
133 const FluxVarCache& fluxVarCache)
134 {
135 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
136 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
137 const auto& insideVolVars = elemVolVars[insideScv];
138 const auto& outsideVolVars = elemVolVars[outsideScv];
139
140 auto insideK = insideVolVars.permeability();
141 auto outsideK = outsideVolVars.permeability();
142
143 // scale with correct extrusion factor
144 insideK *= insideVolVars.extrusionFactor();
145 outsideK *= outsideVolVars.extrusionFactor();
146
147 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
148
149 std::vector<Scalar> ti(fvGeometry.numScv());
150 for (const auto& scv : scvs(fvGeometry))
151 ti[scv.indexInElement()] =
152 -1.0*scvf.area()*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
153
154 return ti;
155 }
156};
157
158} // end namespace Dumux
159
160#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
Darcy's law for box schemes.
Definition: flux/box/darcyslaw.hh:63
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:128
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Definition: flux/box/darcyslaw.hh:79
Declares all properties used in Dumux.