24#ifndef DUMUX_POROUS_MEDIUM_FLOW_FV_SPATIAL_PARAMS_HH
25#define DUMUX_POROUS_MEDIUM_FLOW_FV_SPATIAL_PARAMS_HH
27#include <dune/common/exceptions.hh>
28#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>()))
56template<
class Gr
idGeometry,
class Scalar,
class Implementation>
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;
67 enum { dim = GridView::dimension };
68 enum { dimWorld = GridView::dimensionworld };
69 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
71 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
83 forchCoeffDefault_ = getParam<Scalar>(
"SpatialParams.ForchCoeff", 0.55);
94 [[deprecated(
"Use Dumux::faceTensorAverage from dumux/flux/facetensoraverage.hh instead. This function will be removed after 3.5.")]]
97 const GlobalPosition&
normal)
const
109 [[deprecated(
"Use Dumux::faceTensorAverage from dumux/flux/facetensoraverage.hh instead. This function will be removed after 3.5.")]]
111 const DimWorldMatrix& T2,
112 const GlobalPosition&
normal)
const
126 template<
class ElementSolution>
128 const SubControlVolume& scv,
129 const ElementSolution& elemSol)
const
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");
140 return this->
asImp_().permeabilityAtPos(scv.center());
160 DUNE_THROW(Dune::InvalidStateException,
161 "The spatial parameters do not provide a beaversJosephCoeffAtPos() method.");
172 return forchCoeffDefault_;
176 Scalar forchCoeffDefault_;
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
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