13#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
14#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
16#include <dune/common/reservedvector.hh>
17#include <dune/geometry/type.hh>
18#include <dune/geometry/multilineargeometry.hh>
37template<
class Gr
idView>
40 using Grid =
typename GridView::Grid;
42 static const int dim = Grid::dimension;
43 static const int dimWorld = Grid::dimensionworld;
48 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, GeometryTraits>;
67 using Geometry =
typename T::Geometry;
68 using GridIndexType =
typename T::GridIndexType;
69 using LocalIndexType =
typename T::LocalIndexType;
70 using Scalar =
typename T::Scalar;
71 using GlobalPosition =
typename T::GlobalPosition;
72 using CornerStorage =
typename T::CornerStorage;
73 enum { dim = Geometry::mydimension };
75 static_assert(dim == 2 || dim == 3,
"Box-Dfm sub-control volume only implemented in 2d or 3d");
85 template<
class GeometryHelper>
87 LocalIndexType scvIdx,
90 : isFractureScv_(false)
91 , corners_(geometryHelper.getScvCorners(scvIdx))
94 Dune::GeometryTypes::cube(T::dim),
95 [&](unsigned int i){
return corners_[i]; })
99 , elemLocalScvIdx_(scvIdx)
104 for (
const auto& corner : corners_)
106 center_ /= corners_.size();
118 template<
class GeometryHelper,
class Intersection>
120 const Intersection& intersection,
121 const typename Intersection::Geometry& isGeometry,
122 LocalIndexType indexInIntersection,
123 LocalIndexType vIdxLocal,
124 LocalIndexType elemLocalScvIdx,
125 LocalIndexType elemLocalFacetIdx,
128 : isFractureScv_(true)
133 , vIdxLocal_(vIdxLocal)
134 , elemLocalScvIdx_(elemLocalScvIdx)
136 , facetIdx_(elemLocalFacetIdx)
139 auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
143 corners_[i] = corners[i];
147 Dune::GeometryTypes::cube(T::dim-1),
148 [&](
unsigned int i){
return corners_[i]; }
150 for (
const auto& corner : corners_)
152 center_ /= corners_.size();
165 {
return vIdxLocal_; }
169 {
return elemLocalScvIdx_; }
173 { assert(isFractureScv_);
return facetIdx_; }
177 {
return dofIndex_; }
181 {
return corners_[0]; }
185 {
return elementIndex_; }
189 {
return isFractureScv_; }
193 [[deprecated(
"Will be removed after 3.7. Use fvGeometry.geometry(scv).")]]
197 DUNE_THROW(Dune::InvalidStateException,
"The geometry object cannot be defined for fracture scvs "
198 "because the number of known corners is insufficient. "
199 "You can do this manually by extract the corners from this scv "
200 "and extrude them by the corresponding aperture. ");
202 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
207 CornerStorage corners_;
208 GlobalPosition center_;
210 GridIndexType elementIndex_;
211 LocalIndexType vIdxLocal_;
212 LocalIndexType elemLocalScvIdx_;
213 GridIndexType dofIndex_;
216 LocalIndexType facetIdx_;
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:156
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:172
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:168
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:194
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:188
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:160
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:180
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:79
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:184
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, const Intersection &intersection, const typename Intersection::Geometry &isGeometry, LocalIndexType indexInIntersection, LocalIndexType vIdxLocal, LocalIndexType elemLocalScvIdx, LocalIndexType elemLocalFacetIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Constructor for fracture scvs.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:119
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:164
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:86
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:176
Base class for a sub control volume, i.e a part of the control volume we are making the balance for....
Definition: subcontrolvolumebase.hh:26
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:41
Define some often used mathematical functions.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: common/pdesolver.hh:24
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:39
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:43
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:42
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:40
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:44
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:50
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, GeometryTraits > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:48
typename GeometryTraits::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:49
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:45
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:46
Definition: porousmediumflow/boxdfm/geometryhelper.hh:26
Base class for a sub control volume.
Compute the volume of several common geometry types.