3.5-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>
37
38namespace Dumux {
39
40// forward declaration
41template<class TypeTag, class DiscretizationMethod>
43
44// forward declaration
45template<class Scalar, class GridGeometry>
46class BoxDarcysLaw;
47
52template<class TypeTag>
53class DarcysLawImplementation<TypeTag, DiscretizationMethods::Box>
54: public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
55{ };
56
63template<class Scalar, class GridGeometry>
65{
66 using FVElementGeometry = typename GridGeometry::LocalView;
67 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
68 using Extrusion = Extrusion_t<GridGeometry>;
69 using GridView = typename GridGeometry::GridView;
70 using Element = typename GridView::template Codim<0>::Entity;
71
72 enum { dim = GridView::dimension};
73 enum { dimWorld = GridView::dimensionworld};
74
75public:
76
87 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
88 static Scalar flux(const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace& scvf,
93 const int phaseIdx,
94 const ElementFluxVarsCache& elemFluxVarCache)
95 {
96 const auto& fluxVarCache = elemFluxVarCache[scvf];
97 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
98 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
99 const auto& insideVolVars = elemVolVars[insideScv];
100 const auto& outsideVolVars = elemVolVars[outsideScv];
101
102 auto insideK = insideVolVars.permeability();
103 auto outsideK = outsideVolVars.permeability();
104
105 // scale with correct extrusion factor
106 insideK *= insideVolVars.extrusionFactor();
107 outsideK *= outsideVolVars.extrusionFactor();
108
109 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
110 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
111
112 const auto& shapeValues = fluxVarCache.shapeValues();
113
114 // evaluate gradP - rho*g at integration point
115 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
116 Scalar rho(0.0);
117 for (auto&& scv : scvs(fvGeometry))
118 {
119 const auto& volVars = elemVolVars[scv];
120
121 if (enableGravity)
122 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
123
124 // the global shape function gradient
125 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
126 }
127
128 if (enableGravity)
129 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
130
131 // apply the permeability and return the flux
132 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(scvf);
133 }
134
135 // compute transmissibilities ti for analytical jacobians
136 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
137 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf,
142 const FluxVarCache& fluxVarCache)
143 {
144 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
145 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
146 const auto& insideVolVars = elemVolVars[insideScv];
147 const auto& outsideVolVars = elemVolVars[outsideScv];
148
149 auto insideK = insideVolVars.permeability();
150 auto outsideK = outsideVolVars.permeability();
151
152 // scale with correct extrusion factor
153 insideK *= insideVolVars.extrusionFactor();
154 outsideK *= outsideVolVars.extrusionFactor();
155
156 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
157
158 std::vector<Scalar> ti(fvGeometry.numScv());
159 for (const auto& scv : scvs(fvGeometry))
160 ti[scv.indexInElement()] =
161 -1.0*Extrusion::area(scvf)*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
162
163 return ti;
164 }
165};
166
167} // end namespace Dumux
168
169#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
A free function to average a Tensor at an interface.
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
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:41
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:42
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:65
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:137
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:88
Declares all properties used in Dumux.