25#ifndef DUMUX_SEQUENTIAL_FV_SPATIAL_PARAMS_ONE_P_HH
26#define DUMUX_SEQUENTIAL_FV_SPATIAL_PARAMS_ONE_P_HH
31#include <dune/common/fmatrix.hh>
40template<
class TypeTag>
50 dim = GridView::dimension,
51 dimWorld = GridView::dimensionworld
54 using Element =
typename GridView::template Codim<0>::Entity;
55 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
56 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
67 Scalar
meanK(Scalar K1, Scalar K2)
const
79 void meanK(DimWorldMatrix &result, Scalar K1, Scalar K2)
const
82 for (
int i = 0; i < dimWorld; ++i)
84 for (
int j = 0; j < dimWorld; ++j)
96 void meanK(DimWorldMatrix &result,
const DimWorldMatrix &K1,
const DimWorldMatrix &K2)
const
99 for (
int i = 0; i < dimWorld; ++i)
102 for (
int j = 0; j < dimWorld; ++j)
106 result[i][j] = 0.5 * (K1[i][j] + K2[i][j]);
117 void meanK(DimWorldMatrix &result, Scalar K)
const
119 for (
int i = 0; i < dimWorld; ++i)
121 for (
int j = 0; j < dimWorld; ++j)
132 void meanK(DimWorldMatrix &result,
const DimWorldMatrix &K)
const
145 return asImp_().intrinsicPermeabilityAtPos(element.geometry().center());
156 DUNE_THROW(Dune::InvalidStateException,
157 "The spatial parameters do not provide "
158 "a intrinsicPermeabilityAtPos() method.");
169 return asImp_().porosityAtPos(element.geometry().center());
180 DUNE_THROW(Dune::InvalidStateException,
181 "The spatial parameters do not provide "
182 "a porosityAtPos() method.");
186 Implementation &asImp_()
188 return *
static_cast<Implementation*
> (
this);
191 const Implementation &asImp_()
const
193 return *
static_cast<const Implementation*
> (
this);
Define some often used mathematical functions.
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
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
The base class for spatial parameters of problems using the fv method.
Definition: sequentialfv1p.hh:42
Scalar porosity(const Element &element) const
Function for defining the porosity.
Definition: sequentialfv1p.hh:167
Scalar porosityAtPos(const GlobalPosition &globalPos) const
Function for defining the porosity.
Definition: sequentialfv1p.hh:178
SequentialFVSpatialParamsOneP(const Problem &problem)
Definition: sequentialfv1p.hh:59
const DimWorldMatrix & intrinsicPermeability(const Element &element) const
Function for defining the intrinsic (absolute) permeability.
Definition: sequentialfv1p.hh:143
void meanK(DimWorldMatrix &result, const DimWorldMatrix &K1, const DimWorldMatrix &K2) const
Averages the intrinsic permeability (Tensor).
Definition: sequentialfv1p.hh:96
const DimWorldMatrix & intrinsicPermeabilityAtPos(const GlobalPosition &globalPos) const
Function for defining the intrinsic (absolute) permeability.
Definition: sequentialfv1p.hh:154
void meanK(DimWorldMatrix &result, const DimWorldMatrix &K) const
Dummy function that can be used if only one value exist (boundaries).
Definition: sequentialfv1p.hh:132
Scalar meanK(Scalar K1, Scalar K2) const
Averages the intrinsic permeability (Scalar).
Definition: sequentialfv1p.hh:67
void meanK(DimWorldMatrix &result, Scalar K1, Scalar K2) const
Averages the intrinsic permeability (Scalar).
Definition: sequentialfv1p.hh:79
void meanK(DimWorldMatrix &result, Scalar K) const
Dummy function that can be used if only one value exist (boundaries).
Definition: sequentialfv1p.hh:117
Declares all properties used in Dumux.