3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fv1p.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_FV_SPATIAL_PARAMS_ONE_P_HH
26#define DUMUX_FV_SPATIAL_PARAMS_ONE_P_HH
27
28#include <dune/common/exceptions.hh>
29#include <dune/common/fmatrix.hh>
30
32#include <dumux/common/math.hh>
34
35namespace Dumux {
36
37#ifndef DOXYGEN
38namespace Detail {
39// helper struct detecting if the user-defined spatial params class has a permeabilityAtPos function
40// for g++ > 5.3, this can be replaced by a lambda
41template<class GlobalPosition>
42struct hasPermeabilityAtPos
43{
44 template<class SpatialParams>
45 auto operator()(const SpatialParams& a)
46 -> decltype(a.permeabilityAtPos(std::declval<GlobalPosition>()))
47 {}
48};
49
50template<class GlobalPosition, class SolidSystem>
51struct hasInertVolumeFractionAtPos
52{
53 template<class SpatialParams>
54 auto operator()(const SpatialParams& a)
55 -> decltype(a.template inertVolumeFractionAtPos<SolidSystem>(std::declval<GlobalPosition>(), 0))
56 {}
57};
58
59template<class GlobalPosition>
60struct hasPorosityAtPos
61{
62 template<class SpatialParams>
63 auto operator()(const SpatialParams& a)
64 -> decltype(a.porosityAtPos(std::declval<GlobalPosition>()))
65 {}
66};
67} // end namespace Detail
68#endif
69
75template<class GridGeometry, class Scalar, class Implementation>
77{
78 using GridView = typename GridGeometry::GridView;
79 using FVElementGeometry = typename GridGeometry::LocalView;
80 using SubControlVolume = typename GridGeometry::SubControlVolume;
81 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
82 using Element = typename GridView::template Codim<0>::Entity;
83
84 enum { dim = GridView::dimension };
85 enum { dimWorld = GridView::dimensionworld };
86 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
87
88 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
89
90public:
91 FVSpatialParamsOneP(std::shared_ptr<const GridGeometry> gridGeometry)
92 : gridGeometry_(gridGeometry)
93 , gravity_(0.0)
94 {
95 const bool enableGravity = getParam<bool>("Problem.EnableGravity");
96 if (enableGravity)
97 gravity_[dimWorld-1] = -9.81;
98
99 /* \brief default forchheimer coefficient
100 * Source: Ward, J.C. 1964 Turbulent flow in porous media. ASCE J. Hydraul. Div 90 \cite ward1964 .
101 * Actually the Forchheimer coefficient is also a function of the dimensions of the
102 * porous medium. Taking it as a constant is only a first approximation
103 * (Nield, Bejan, Convection in porous media, 2006, p. 10 \cite nield2006 )
104 */
105 forchCoeffDefault_ = getParam<Scalar>("SpatialParams.ForchCoeff", 0.55);
106 }
107
118 const GlobalPosition& gravity(const GlobalPosition &pos) const
119 { return gravity_; }
120
129 Scalar harmonicMean(const Scalar T1,
130 const Scalar T2,
131 const GlobalPosition& normal) const
132 { return Dumux::harmonicMean(T1, T2); }
133
143 DimWorldMatrix harmonicMean(const DimWorldMatrix& T1,
144 const DimWorldMatrix& T2,
145 const GlobalPosition& normal) const
146 {
147 // determine nT*k*n
148 GlobalPosition tmp;
149 GlobalPosition tmp2;
150 T1.mv(normal, tmp);
151 T2.mv(normal, tmp2);
152 const Scalar alpha1 = tmp*normal;
153 const Scalar alpha2 = tmp2*normal;
154
155 const Scalar alphaHarmonic = Dumux::harmonicMean(alpha1, alpha2);
156 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
157
158 DimWorldMatrix T(0.0);
159 for (int i = 0; i < dimWorld; ++i)
160 {
161 for (int j = 0; j < dimWorld; ++j)
162 {
163 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
164 if (i == j)
165 T[i][j] += alphaHarmonic - alphaAverage;
166 }
167 }
168
169 return T;
170 }
171
181 template<class ElementSolution>
182 decltype(auto) permeability(const Element& element,
183 const SubControlVolume& scv,
184 const ElementSolution& elemSol) const
185 {
186 static_assert(decltype(isValid(Detail::hasPermeabilityAtPos<GlobalPosition>())(this->asImp_()))::value," \n\n"
187 " Your spatial params class has to either implement\n\n"
188 " const PermeabilityType& permeabilityAtPos(const GlobalPosition& globalPos) const\n\n"
189 " or overload this function\n\n"
190 " template<class ElementSolution>\n"
191 " const PermeabilityType& permeability(const Element& element,\n"
192 " const SubControlVolume& scv,\n"
193 " const ElementSolution& elemSol) const\n\n");
194
195 return asImp_().permeabilityAtPos(scv.center());
196 }
197
203 static constexpr bool evaluatePermeabilityAtScvfIP()
204 { return false; }
205
216 template<class ElementSolution>
217 Scalar porosity(const Element& element,
218 const SubControlVolume& scv,
219 const ElementSolution& elemSol) const
220 {
221 static_assert(decltype(isValid(Detail::hasPorosityAtPos<GlobalPosition>())(this->asImp_()))::value," \n\n"
222 " Your spatial params class has to either implement\n\n"
223 " Scalar porosityAtPos(const GlobalPosition& globalPos) const\n\n"
224 " or overload this function\n\n"
225 " template<class ElementSolution>\n"
226 " Scalar porosity(const Element& element,\n"
227 " const SubControlVolume& scv,\n"
228 " const ElementSolution& elemSol) const\n\n");
229
230 return asImp_().porosityAtPos(scv.center());
231 }
232
249 template<class SolidSystem, class ElementSolution,
250 typename std::enable_if_t<SolidSystem::isInert()
251 && SolidSystem::numInertComponents == 1
252 && !decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())(std::declval<Implementation>()))::value,
253 int> = 0>
254 Scalar inertVolumeFraction(const Element& element,
255 const SubControlVolume& scv,
256 const ElementSolution& elemSol,
257 int compIdx) const
258 {
259 return 1.0 - asImp_().porosity(element, scv, elemSol);
260 }
261
262 // specialization if there are no inert components at all
263 template<class SolidSystem, class ElementSolution,
264 typename std::enable_if_t<SolidSystem::numInertComponents == 0, int> = 0>
265 Scalar inertVolumeFraction(const Element& element,
266 const SubControlVolume& scv,
267 const ElementSolution& elemSol,
268 int compIdx) const
269 {
270 return 0.0;
271 }
272
273 // the more general interface forwarding to inertVolumeFractionAtPos
274 template<class SolidSystem, class ElementSolution,
275 typename std::enable_if_t<(SolidSystem::numInertComponents > 1) ||
276 (
277 (SolidSystem::numInertComponents > 0) &&
278 (
279 !SolidSystem::isInert()
280 || decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())
281 (std::declval<Implementation>()))::value
282 )
283 ),
284 int> = 0>
285 Scalar inertVolumeFraction(const Element& element,
286 const SubControlVolume& scv,
287 const ElementSolution& elemSol,
288 int compIdx) const
289 {
290 static_assert(decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())(this->asImp_()))::value," \n\n"
291 " Your spatial params class has to either implement\n\n"
292 " template<class SolidSystem>\n"
293 " Scalar inertVolumeFractionAtPos(const GlobalPosition& globalPos, int compIdx) const\n\n"
294 " or overload this function\n\n"
295 " template<class SolidSystem, class ElementSolution>\n"
296 " Scalar inertVolumeFraction(const Element& element,\n"
297 " const SubControlVolume& scv,\n"
298 " const ElementSolution& elemSol,\n"
299 " int compIdx) const\n\n");
300
301 return asImp_().template inertVolumeFractionAtPos<SolidSystem>(scv.center(), compIdx);
302 }
303
311 Scalar beaversJosephCoeffAtPos(const GlobalPosition& globalPos) const
312 {
313 DUNE_THROW(Dune::InvalidStateException,
314 "The spatial parameters do not provide a beaversJosephCoeffAtPos() method.");
315 }
316
323 Scalar forchCoeff(const SubControlVolumeFace &scvf) const
324 {
325 return forchCoeffDefault_;
326 }
327
329 const GridGeometry& gridGeometry() const
330 { return *gridGeometry_; }
331
332
333protected:
334 Implementation &asImp_()
335 { return *static_cast<Implementation*>(this); }
336
337 const Implementation &asImp_() const
338 { return *static_cast<const Implementation*>(this); }
339
340private:
341 std::shared_ptr<const GridGeometry> gridGeometry_;
342 GlobalPosition gravity_;
343 Scalar forchCoeffDefault_;
344};
345
346} // namespace Dumux
347
348#endif
A helper function for class member function introspection.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
Definition: adapt.hh:29
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
The base class for spatial parameters of one-phase problems using a fully implicit discretization met...
Definition: fv1p.hh:77
DimWorldMatrix harmonicMean(const DimWorldMatrix &T1, const DimWorldMatrix &T2, const GlobalPosition &normal) const
Harmonic average of a discontinuous tensorial field at discontinuity interface.
Definition: fv1p.hh:143
Scalar harmonicMean(const Scalar T1, const Scalar T2, const GlobalPosition &normal) const
Harmonic average of a discontinuous scalar field at discontinuity interface (for compatibility reason...
Definition: fv1p.hh:129
Scalar beaversJosephCoeffAtPos(const GlobalPosition &globalPos) const
Function for defining the Beavers-Joseph coefficient for multidomain problems .
Definition: fv1p.hh:311
static constexpr bool evaluatePermeabilityAtScvfIP()
If the permeability should be evaluated directly at the scvf integration point (for convergence tests...
Definition: fv1p.hh:203
const Implementation & asImp_() const
Definition: fv1p.hh:337
Implementation & asImp_()
Definition: fv1p.hh:334
Scalar forchCoeff(const SubControlVolumeFace &scvf) const
Apply the Forchheimer coefficient for inertial forces calculation.
Definition: fv1p.hh:323
const GlobalPosition & gravity(const GlobalPosition &pos) const
Returns the acceleration due to gravity .
Definition: fv1p.hh:118
Scalar porosity(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Function for defining the porosity. That is possibly solution dependent.
Definition: fv1p.hh:217
FVSpatialParamsOneP(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: fv1p.hh:91
Scalar inertVolumeFraction(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol, int compIdx) const
Function for defining the solid volume fraction. That is possibly solution dependent.
Definition: fv1p.hh:254
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: fv1p.hh:329
decltype(auto) permeability(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Function for defining the (intrinsic) permeability .
Definition: fv1p.hh:182