25#ifndef DUMUX_FV_SPATIAL_PARAMS_ONE_P_HH
26#define DUMUX_FV_SPATIAL_PARAMS_ONE_P_HH
28#include <dune/common/exceptions.hh>
29#include <dune/common/fmatrix.hh>
41template<
class GlobalPosition>
42struct hasPermeabilityAtPos
44 template<
class SpatialParams>
45 auto operator()(
const SpatialParams& a)
46 ->
decltype(a.permeabilityAtPos(std::declval<GlobalPosition>()))
50template<
class GlobalPosition,
class Sol
idSystem>
51struct hasInertVolumeFractionAtPos
53 template<
class SpatialParams>
54 auto operator()(
const SpatialParams& a)
55 ->
decltype(a.template inertVolumeFractionAtPos<SolidSystem>(std::declval<GlobalPosition>(), 0))
59template<
class GlobalPosition>
60struct hasPorosityAtPos
62 template<
class SpatialParams>
63 auto operator()(
const SpatialParams& a)
64 ->
decltype(a.porosityAtPos(std::declval<GlobalPosition>()))
75template<
class Gr
idGeometry,
class Scalar,
class Implementation>
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;
84 enum { dim = GridView::dimension };
85 enum { dimWorld = GridView::dimensionworld };
86 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
88 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
95 const bool enableGravity = getParam<bool>(
"Problem.EnableGravity");
97 gravity_[dimWorld-1] = -9.81;
105 forchCoeffDefault_ = getParam<Scalar>(
"SpatialParams.ForchCoeff", 0.55);
118 const GlobalPosition&
gravity(
const GlobalPosition &pos)
const
131 const GlobalPosition& normal)
const
144 const DimWorldMatrix& T2,
145 const GlobalPosition& normal)
const
152 const Scalar alpha1 = tmp*normal;
153 const Scalar alpha2 = tmp2*normal;
156 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
158 DimWorldMatrix T(0.0);
159 for (
int i = 0; i < dimWorld; ++i)
161 for (
int j = 0; j < dimWorld; ++j)
163 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
165 T[i][j] += alphaHarmonic - alphaAverage;
181 template<
class ElementSolution>
183 const SubControlVolume& scv,
184 const ElementSolution& elemSol)
const
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");
195 return asImp_().permeabilityAtPos(scv.center());
216 template<
class ElementSolution>
218 const SubControlVolume& scv,
219 const ElementSolution& elemSol)
const
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");
230 return asImp_().porosityAtPos(scv.center());
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,
255 const SubControlVolume& scv,
256 const ElementSolution& elemSol,
259 return 1.0 -
asImp_().porosity(element, scv, elemSol);
263 template<
class SolidSystem,
class ElementSolution,
264 typename std::enable_if_t<SolidSystem::numInertComponents == 0, int> = 0>
266 const SubControlVolume& scv,
267 const ElementSolution& elemSol,
274 template<
class SolidSystem,
class ElementSolution,
275 typename std::enable_if_t<(SolidSystem::numInertComponents > 1) ||
277 (SolidSystem::numInertComponents > 0) &&
279 !SolidSystem::isInert()
280 ||
decltype(
isValid(Detail::hasInertVolumeFractionAtPos<GlobalPosition, SolidSystem>())
281 (std::declval<Implementation>()))::value
286 const SubControlVolume& scv,
287 const ElementSolution& elemSol,
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");
301 return asImp_().template inertVolumeFractionAtPos<SolidSystem>(scv.center(), compIdx);
313 DUNE_THROW(Dune::InvalidStateException,
314 "The spatial parameters do not provide a beaversJosephCoeffAtPos() method.");
325 return forchCoeffDefault_;
330 {
return *gridGeometry_; }
335 {
return *
static_cast<Implementation*
>(
this); }
338 {
return *
static_cast<const Implementation*
>(
this); }
341 std::shared_ptr<const GridGeometry> gridGeometry_;
342 GlobalPosition gravity_;
343 Scalar forchCoeffDefault_;
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
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