25#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
28#include <dune/common/reservedvector.hh>
29#include <dune/geometry/type.hh>
30#include <dune/geometry/multilineargeometry.hh>
49template<
class Gr
idView>
52 using Grid =
typename GridView::Grid;
54 static const int dim = Grid::dimension;
55 static const int dimWorld = Grid::dimensionworld;
60 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, GeometryTraits>;
79 using Geometry =
typename T::Geometry;
80 using GridIndexType =
typename T::GridIndexType;
81 using LocalIndexType =
typename T::LocalIndexType;
82 using Scalar =
typename T::Scalar;
83 using GlobalPosition =
typename T::GlobalPosition;
84 using CornerStorage =
typename T::CornerStorage;
85 enum { dim = Geometry::mydimension };
87 static_assert(dim == 2 || dim == 3,
"Box-Dfm sub-control volume only implemented in 2d or 3d");
97 template<
class GeometryHelper>
99 LocalIndexType scvIdx,
102 : isFractureScv_(false)
103 , corners_(geometryHelper.getScvCorners(scvIdx))
106 Dune::GeometryTypes::cube(T::dim),
107 [&](unsigned int i){
return corners_[i]; })
111 , elemLocalScvIdx_(scvIdx)
116 for (
const auto&
corner : corners_)
118 center_ /= corners_.size();
130 template<
class GeometryHelper,
class Intersection>
132 const Intersection& intersection,
133 const typename Intersection::Geometry& isGeometry,
134 LocalIndexType indexInIntersection,
135 LocalIndexType vIdxLocal,
136 LocalIndexType elemLocalScvIdx,
137 LocalIndexType elemLocalFacetIdx,
140 : isFractureScv_(true)
145 , vIdxLocal_(vIdxLocal)
146 , elemLocalScvIdx_(elemLocalScvIdx)
148 , facetIdx_(elemLocalFacetIdx)
151 auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
155 corners_[i] = corners[i];
159 Dune::GeometryTypes::cube(T::dim-1),
160 [&](
unsigned int i){
return corners_[i]; }
162 for (
const auto&
corner : corners_)
164 center_ /= corners_.size();
177 {
return vIdxLocal_; }
181 {
return elemLocalScvIdx_; }
185 { assert(isFractureScv_);
return facetIdx_; }
189 {
return dofIndex_; }
193 {
return corners_[0]; }
197 {
return elementIndex_; }
201 {
return isFractureScv_; }
208 DUNE_THROW(Dune::InvalidStateException,
"The geometry object cannot be defined for fracture scvs "
209 "because the number of known corners is insufficient. "
210 "You can do this manually by extract the corners from this scv "
211 "and extrude them by the corresponding aperture. ");
213 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
217 const GlobalPosition&
corner(LocalIndexType localIdx)
const
219 assert(localIdx < corners_.size() &&
"provided index exceeds the number of corners");
220 return corners_[localIdx];
225 CornerStorage corners_;
226 GlobalPosition center_;
228 GridIndexType elementIndex_;
229 LocalIndexType vIdxLocal_;
230 LocalIndexType elemLocalScvIdx_;
231 GridIndexType dofIndex_;
234 LocalIndexType facetIdx_;
Define some often used mathematical functions.
Base class for a sub control volume.
Compute the volume of several common geometry types.
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:53
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
Base class for a sub control volume, i.e a part of the control volume we are making the balance for....
Definition: subcontrolvolumebase.hh:38
Definition: porousmediumflow/boxdfm/geometryhelper.hh:38
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:51
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:55
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:54
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:52
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:56
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:62
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, GeometryTraits > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:60
typename GeometryTraits::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:61
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:57
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:58
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:76
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:168
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:184
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:180
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:205
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:200
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:217
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:172
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:192
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:91
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:196
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:131
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:176
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:98
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:188
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.