25#ifndef DUMUX_INCOMPRESSIBLE_ONEP_TEST_SPATIAL_PARAMS_HH
26#define DUMUX_INCOMPRESSIBLE_ONEP_TEST_SPATIAL_PARAMS_HH
32namespace LensSpatialParams {
40template<
class GlobalPosition>
42 const GlobalPosition& lensLowerLeft,
43 const GlobalPosition& lensUpperRight)
45 const auto eps = 1e-8*(lensUpperRight - lensLowerLeft).two_norm();
46 for (
int i = 0; i < GlobalPosition::size(); ++i)
47 if (globalPos[i] < lensLowerLeft[i] +
eps || globalPos[i] > lensUpperRight[i] -
eps)
59template<
class Gr
idGeometry,
class Scalar>
61:
public FVSpatialParamsOneP<GridGeometry, Scalar, OnePTestSpatialParams<GridGeometry, Scalar>>
65 using GridView =
typename GridGeometry::GridView;
66 using Element =
typename GridView::template Codim<0>::Entity;
67 using FVElementGeometry =
typename GridGeometry::LocalView;
68 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
70 static constexpr int dimWorld = GridView::dimensionworld;
71 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
77 , lensLowerLeft_(std::numeric_limits<Scalar>::max())
78 , lensUpperRight_(std::numeric_limits<Scalar>::lowest())
80 permeability_ = getParam<Scalar>(
"SpatialParams.Permeability");
81 permeabilityLens_ = getParam<Scalar>(
"SpatialParams.PermeabilityLens");
90 {
return isInLens_(globalPos) ? permeabilityLens_ : permeability_; }
103 void setLens(
const GlobalPosition& lowerLeft,
const GlobalPosition& upperRight)
104 { lensLowerLeft_ = lowerLeft; lensUpperRight_ = upperRight; }
107 bool isInLens_(
const GlobalPosition &globalPos)
const
110 GlobalPosition lensLowerLeft_, lensUpperRight_;
111 Scalar permeability_, permeabilityLens_;
The base class for spatial parameters of one-phase problems using a fully implicit discretization met...
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
bool pointInLens(const GlobalPosition &globalPos, const GlobalPosition &lensLowerLeft, const GlobalPosition &lensUpperRight)
Returns if a point is in a lens with a given bounding box.
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:41
The base class for spatial parameters of one-phase problems using a fully implicit discretization met...
Definition: fv1p.hh:77
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: fv1p.hh:334
The spatial parameters class for the test problem using the incompressible 1p model.
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:62
Scalar PermeabilityType
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:74
void setLens(const GlobalPosition &lowerLeft, const GlobalPosition &upperRight)
Optionally set a lens.
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:103
PermeabilityType permeabilityAtPos(const GlobalPosition &globalPos) const
Function for defining the (intrinsic) permeability .
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:89
Scalar porosityAtPos(const GlobalPosition &globalPos) const
Defines the porosity .
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:97
OnePTestSpatialParams(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:75