3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/fvspatialparams.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 *****************************************************************************/
24#ifndef DUMUX_POROUS_MEDIUM_FLOW_FV_SPATIAL_PARAMS_HH
25#define DUMUX_POROUS_MEDIUM_FLOW_FV_SPATIAL_PARAMS_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/common/fmatrix.hh>
29
31#include <dumux/common/math.hh>
32#include <dumux/flux/facetensoraverage.hh> // remove include after release 3.5
35
36namespace Dumux {
37
38#ifndef DOXYGEN
39namespace Detail {
40// helper struct detecting if the user-defined spatial params class has a permeabilityAtPos function
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} // end namespace Detail
50#endif
51
56template<class GridGeometry, class Scalar, class Implementation>
58: public FVPorousMediumSpatialParams<GridGeometry, Scalar, Implementation>
59{
61 using GridView = typename GridGeometry::GridView;
62 using FVElementGeometry = typename GridGeometry::LocalView;
63 using SubControlVolume = typename GridGeometry::SubControlVolume;
64 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
65 using Element = typename GridView::template Codim<0>::Entity;
66
67 enum { dim = GridView::dimension };
68 enum { dimWorld = GridView::dimensionworld };
69 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
70
71 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
72
73public:
74 FVPorousMediumFlowSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry)
76 {
77 /* \brief default forchheimer coefficient
78 * Source: Ward, J.C. 1964 Turbulent flow in porous media. ASCE J. Hydraul. Div 90 \cite ward1964 .
79 * Actually the Forchheimer coefficient is also a function of the dimensions of the
80 * porous medium. Taking it as a constant is only a first approximation
81 * (Nield, Bejan, Convection in porous media, 2006, p. 10 \cite nield2006 )
82 */
83 forchCoeffDefault_ = getParam<Scalar>("SpatialParams.ForchCoeff", 0.55);
84 }
85
94 [[deprecated("Use Dumux::faceTensorAverage from dumux/flux/facetensoraverage.hh instead. This function will be removed after 3.5.")]]
95 Scalar harmonicMean(const Scalar T1,
96 const Scalar T2,
97 const GlobalPosition& normal) const
98 { return Dumux::harmonicMean(T1, T2); }
99
109 [[deprecated("Use Dumux::faceTensorAverage from dumux/flux/facetensoraverage.hh instead. This function will be removed after 3.5.")]]
110 DimWorldMatrix harmonicMean(const DimWorldMatrix& T1,
111 const DimWorldMatrix& T2,
112 const GlobalPosition& normal) const
113 {
114 return faceTensorAverage(T1, T2, normal);
115 }
116
126 template<class ElementSolution>
127 decltype(auto) permeability(const Element& element,
128 const SubControlVolume& scv,
129 const ElementSolution& elemSol) const
130 {
131 static_assert(decltype(isValid(Detail::hasPermeabilityAtPos<GlobalPosition>())(this->asImp_()))::value," \n\n"
132 " Your spatial params class has to either implement\n\n"
133 " const PermeabilityType& permeabilityAtPos(const GlobalPosition& globalPos) const\n\n"
134 " or overload this function\n\n"
135 " template<class ElementSolution>\n"
136 " const PermeabilityType& permeability(const Element& element,\n"
137 " const SubControlVolume& scv,\n"
138 " const ElementSolution& elemSol) const\n\n");
139
140 return this->asImp_().permeabilityAtPos(scv.center());
141 }
142
148 static constexpr bool evaluatePermeabilityAtScvfIP()
149 { return false; }
150
158 Scalar beaversJosephCoeffAtPos(const GlobalPosition& globalPos) const
159 {
160 DUNE_THROW(Dune::InvalidStateException,
161 "The spatial parameters do not provide a beaversJosephCoeffAtPos() method.");
162 }
163
170 Scalar forchCoeff(const SubControlVolumeFace &scvf) const
171 {
172 return forchCoeffDefault_;
173 }
174
175private:
176 Scalar forchCoeffDefault_;
177};
178
179} // namespace Dumux
180
181#endif
A helper function for class member function introspection.
The base class for spatial parameters in porous medium problems.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A free function to average a Tensor at an interface.
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
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
The base class for spatial parameters of porous-medium problems.
Definition: fvporousmediumspatialparams.hh:63
The base class for spatial parameters used with finite-volume schemes.
Definition: common/fvspatialparams.hh:46
Implementation & asImp_()
Returns the implementation of the spatial parameters (static polymorphism)
Definition: common/fvspatialparams.hh:147
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvspatialparams.hh:142
The base class for spatial parameters of porous-medium-flow problems.
Definition: porousmediumflow/fvspatialparams.hh:59
decltype(auto) permeability(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Function for defining the (intrinsic) permeability .
Definition: porousmediumflow/fvspatialparams.hh:127
static constexpr bool evaluatePermeabilityAtScvfIP()
If the permeability should be evaluated directly at the scvf integration point (for convergence tests...
Definition: porousmediumflow/fvspatialparams.hh:148
FVPorousMediumFlowSpatialParams(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: porousmediumflow/fvspatialparams.hh:74
DimWorldMatrix harmonicMean(const DimWorldMatrix &T1, const DimWorldMatrix &T2, const GlobalPosition &normal) const
Harmonic average of a discontinuous tensorial field at discontinuity interface.
Definition: porousmediumflow/fvspatialparams.hh:110
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: porousmediumflow/fvspatialparams.hh:95
Scalar beaversJosephCoeffAtPos(const GlobalPosition &globalPos) const
Function for defining the Beavers-Joseph coefficient for multidomain problems .
Definition: porousmediumflow/fvspatialparams.hh:158
Scalar forchCoeff(const SubControlVolumeFace &scvf) const
Apply the Forchheimer coefficient for inertial forces calculation.
Definition: porousmediumflow/fvspatialparams.hh:170