3.5-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#warning "This file is deprecated, use FVPorousMediumFlowSpatialParamsOneP from dumux/porousmediumflow/fvspatialparams1p.hh instead!"
29
30#include <dune/common/exceptions.hh>
31#include <dune/common/fmatrix.hh>
32
34#include <dumux/common/math.hh>
36
37namespace Dumux {
38
39#ifndef DOXYGEN
40namespace Detail {
41// helper struct detecting if the user-defined spatial params class has a permeabilityAtPos function
42template<class GlobalPosition>
43struct hasPermeabilityAtPos
44{
45 template<class SpatialParams>
46 auto operator()(const SpatialParams& a)
47 -> decltype(a.permeabilityAtPos(std::declval<GlobalPosition>()))
48 {}
49};
50
51template<class GlobalPosition, class SolidSystem>
52struct hasInertVolumeFractionAtPos
53{
54 template<class SpatialParams>
55 auto operator()(const SpatialParams& a)
56 -> decltype(a.template inertVolumeFractionAtPos<SolidSystem>(std::declval<GlobalPosition>(), 0))
57 {}
58};
59
60template<class GlobalPosition>
61struct hasPorosityAtPos
62{
63 template<class SpatialParams>
64 auto operator()(const SpatialParams& a)
65 -> decltype(a.porosityAtPos(std::declval<GlobalPosition>()))
66 {}
67};
68} // end namespace Detail
69#endif
70
76template<class GridGeometry, class Scalar, class Implementation>
77class
78[[deprecated("Use FVPorousMediumFlowSpatialParamsOneP from dumux/porousmediumflow/fvspatialparamsmp.hh instead. This class will be removed after 3.5.")]]
80{
81 using GridView = typename GridGeometry::GridView;
82 using FVElementGeometry = typename GridGeometry::LocalView;
83 using SubControlVolume = typename GridGeometry::SubControlVolume;
84 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
85 using Element = typename GridView::template Codim<0>::Entity;
86
87 enum { dim = GridView::dimension };
88 enum { dimWorld = GridView::dimensionworld };
89 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
90
91 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
92
93public:
94 FVSpatialParamsOneP(std::shared_ptr<const GridGeometry> gridGeometry)
95 : gridGeometry_(gridGeometry)
96 , gravity_(0.0)
97 {
98 const bool enableGravity = getParam<bool>("Problem.EnableGravity");
99 if (enableGravity)
100 gravity_[dimWorld-1] = -9.81;
101
102 /* \brief default forchheimer coefficient
103 * Source: Ward, J.C. 1964 Turbulent flow in porous media. ASCE J. Hydraul. Div 90 \cite ward1964 .
104 * Actually the Forchheimer coefficient is also a function of the dimensions of the
105 * porous medium. Taking it as a constant is only a first approximation
106 * (Nield, Bejan, Convection in porous media, 2006, p. 10 \cite nield2006 )
107 */
108 forchCoeffDefault_ = getParam<Scalar>("SpatialParams.ForchCoeff", 0.55);
109 }
110
121 const GlobalPosition& gravity(const GlobalPosition &pos) const
122 { return gravity_; }
123
132 Scalar harmonicMean(const Scalar T1,
133 const Scalar T2,
134 const GlobalPosition& normal) const
135 { return Dumux::harmonicMean(T1, T2); }
136
146 DimWorldMatrix harmonicMean(const DimWorldMatrix& T1,
147 const DimWorldMatrix& T2,
148 const GlobalPosition& normal) const
149 {
150 // determine nT*k*n
151 GlobalPosition tmp;
152 GlobalPosition tmp2;
153 T1.mv(normal, tmp);
154 T2.mv(normal, tmp2);
155 const Scalar alpha1 = tmp*normal;
156 const Scalar alpha2 = tmp2*normal;
157
158 const Scalar alphaHarmonic = Dumux::harmonicMean(alpha1, alpha2);
159 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
160
161 DimWorldMatrix T(0.0);
162 for (int i = 0; i < dimWorld; ++i)
163 {
164 for (int j = 0; j < dimWorld; ++j)
165 {
166 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
167 if (i == j)
168 T[i][j] += alphaHarmonic - alphaAverage;
169 }
170 }
171
172 return T;
173 }
174
184 template<class ElementSolution>
185 decltype(auto) permeability(const Element& element,
186 const SubControlVolume& scv,
187 const ElementSolution& elemSol) const
188 {
189 static_assert(decltype(isValid(Detail::hasPermeabilityAtPos<GlobalPosition>())(this->asImp_()))::value," \n\n"
190 " Your spatial params class has to either implement\n\n"
191 " const PermeabilityType& permeabilityAtPos(const GlobalPosition& globalPos) const\n\n"
192 " or overload this function\n\n"
193 " template<class ElementSolution>\n"
194 " const PermeabilityType& permeability(const Element& element,\n"
195 " const SubControlVolume& scv,\n"
196 " const ElementSolution& elemSol) const\n\n");
197
198 return asImp_().permeabilityAtPos(scv.center());
199 }
200
206 static constexpr bool evaluatePermeabilityAtScvfIP()
207 { return false; }
208
219 template<class ElementSolution>
220 Scalar porosity(const Element& element,
221 const SubControlVolume& scv,
222 const ElementSolution& elemSol) const
223 {
224 static_assert(decltype(isValid(Detail::hasPorosityAtPos<GlobalPosition>())(this->asImp_()))::value," \n\n"
225 " Your spatial params class has to either implement\n\n"
226 " Scalar porosityAtPos(const GlobalPosition& globalPos) const\n\n"
227 " or overload this function\n\n"
228 " template<class ElementSolution>\n"
229 " Scalar porosity(const Element& element,\n"
230 " const SubControlVolume& scv,\n"
231 " const ElementSolution& elemSol) const\n\n");
232
233 return asImp_().porosityAtPos(scv.center());
234 }
235
252 template<class SolidSystem, class ElementSolution,
253 typename std::enable_if_t<SolidSystem::isInert()
254 && SolidSystem::numInertComponents == 1
255 && !decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())(std::declval<Implementation>()))::value,
256 int> = 0>
257 Scalar inertVolumeFraction(const Element& element,
258 const SubControlVolume& scv,
259 const ElementSolution& elemSol,
260 int compIdx) const
261 {
262 return 1.0 - asImp_().porosity(element, scv, elemSol);
263 }
264
265 // specialization if there are no inert components at all
266 template<class SolidSystem, class ElementSolution,
267 typename std::enable_if_t<SolidSystem::numInertComponents == 0, int> = 0>
268 Scalar inertVolumeFraction(const Element& element,
269 const SubControlVolume& scv,
270 const ElementSolution& elemSol,
271 int compIdx) const
272 {
273 return 0.0;
274 }
275
276 // the more general interface forwarding to inertVolumeFractionAtPos
277 template<class SolidSystem, class ElementSolution,
278 typename std::enable_if_t<(SolidSystem::numInertComponents > 1) ||
279 (
280 (SolidSystem::numInertComponents > 0) &&
281 (
282 !SolidSystem::isInert()
283 || decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())
284 (std::declval<Implementation>()))::value
285 )
286 ),
287 int> = 0>
288 Scalar inertVolumeFraction(const Element& element,
289 const SubControlVolume& scv,
290 const ElementSolution& elemSol,
291 int compIdx) const
292 {
293 static_assert(decltype(isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())(this->asImp_()))::value," \n\n"
294 " Your spatial params class has to either implement\n\n"
295 " template<class SolidSystem>\n"
296 " Scalar inertVolumeFractionAtPos(const GlobalPosition& globalPos, int compIdx) const\n\n"
297 " or overload this function\n\n"
298 " template<class SolidSystem, class ElementSolution>\n"
299 " Scalar inertVolumeFraction(const Element& element,\n"
300 " const SubControlVolume& scv,\n"
301 " const ElementSolution& elemSol,\n"
302 " int compIdx) const\n\n");
303
304 return asImp_().template inertVolumeFractionAtPos<SolidSystem>(scv.center(), compIdx);
305 }
306
314 Scalar beaversJosephCoeffAtPos(const GlobalPosition& globalPos) const
315 {
316 DUNE_THROW(Dune::InvalidStateException,
317 "The spatial parameters do not provide a beaversJosephCoeffAtPos() method.");
318 }
319
326 Scalar forchCoeff(const SubControlVolumeFace &scvf) const
327 {
328 return forchCoeffDefault_;
329 }
330
332 const GridGeometry& gridGeometry() const
333 { return *gridGeometry_; }
334
335
336protected:
337 Implementation &asImp_()
338 { return *static_cast<Implementation*>(this); }
339
340 const Implementation &asImp_() const
341 { return *static_cast<const Implementation*>(this); }
342
343private:
344 std::shared_ptr<const GridGeometry> gridGeometry_;
345 GlobalPosition gravity_;
346 Scalar forchCoeffDefault_;
347};
348
349} // namespace Dumux
350
351#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.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
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:69
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:80
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:146
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:132
Scalar beaversJosephCoeffAtPos(const GlobalPosition &globalPos) const
Function for defining the Beavers-Joseph coefficient for multidomain problems .
Definition: fv1p.hh:314
static constexpr bool evaluatePermeabilityAtScvfIP()
If the permeability should be evaluated directly at the scvf integration point (for convergence tests...
Definition: fv1p.hh:206
const Implementation & asImp_() const
Definition: fv1p.hh:340
Implementation & asImp_()
Definition: fv1p.hh:337
Scalar forchCoeff(const SubControlVolumeFace &scvf) const
Apply the Forchheimer coefficient for inertial forces calculation.
Definition: fv1p.hh:326
const GlobalPosition & gravity(const GlobalPosition &pos) const
Returns the acceleration due to gravity .
Definition: fv1p.hh:121
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:220
FVSpatialParamsOneP(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: fv1p.hh:94
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:257
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: fv1p.hh:332
decltype(auto) permeability(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Function for defining the (intrinsic) permeability .
Definition: fv1p.hh:185