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>
48template<
class Gr
idView>
51 using Grid =
typename GridView::Grid;
53 static const int dim = Grid::dimension;
54 static const int dimWorld = Grid::dimensionworld;
62 template<
int mydim,
int cdim >
65 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(
dim)) >;
72 static const bool v =
true;
73 static const unsigned int topologyId = Dune::GeometryTypes::cube(mydim).id();
80 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, ScvMLGTraits<Scalar>>;
99 using Geometry =
typename T::Geometry;
100 using GridIndexType =
typename T::GridIndexType;
101 using LocalIndexType =
typename T::LocalIndexType;
102 using Scalar =
typename T::Scalar;
103 using GlobalPosition =
typename T::GlobalPosition;
104 using CornerStorage =
typename T::CornerStorage;
105 enum { dim = Geometry::mydimension };
107 static_assert(dim == 2 || dim == 3,
"Box-Dfm sub-control volume only implemented in 2d or 3d");
117 template<
class GeometryHelper>
119 LocalIndexType scvIdx,
122 : isFractureScv_(false)
123 , corners_(geometryHelper.getScvCorners(scvIdx))
125 , volume_(geometryHelper.scvVolume(corners_))
128 , elemLocalScvIdx_(scvIdx)
133 for (
const auto&
corner : corners_)
135 center_ /= corners_.size();
147 template<
class GeometryHelper,
class Intersection>
149 const Intersection& intersection,
150 const typename Intersection::Geometry& isGeometry,
151 LocalIndexType indexInIntersection,
152 LocalIndexType vIdxLocal,
153 LocalIndexType elemLocalScvIdx,
154 LocalIndexType elemLocalFacetIdx,
157 : isFractureScv_(true)
162 , vIdxLocal_(vIdxLocal)
163 , elemLocalScvIdx_(elemLocalScvIdx)
165 , facetIdx_(elemLocalFacetIdx)
168 auto corners = geometryHelper.getBoundaryScvfCorners(intersection, isGeometry, indexInIntersection);
172 corners_[i] = corners[i];
175 volume_ = geometryHelper.scvfArea(corners);
176 for (
const auto&
corner : corners_)
178 center_ /= corners_.size();
191 {
return vIdxLocal_; }
195 {
return elemLocalScvIdx_; }
199 { assert(isFractureScv_);
return facetIdx_; }
203 {
return dofIndex_; }
207 {
return corners_[0]; }
211 {
return elementIndex_; }
215 {
return isFractureScv_; }
222 DUNE_THROW(Dune::InvalidStateException,
"The geometry object cannot be defined for fracture scvs "
223 "because the number of known corners is insufficient. "
224 "You can do this manually by extract the corners from this scv "
225 "and extrude them by the corresponding aperture. ");
227 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
231 const GlobalPosition&
corner(LocalIndexType localIdx)
const
233 assert(localIdx < corners_.size() &&
"provided index exceeds the number of corners");
234 return corners_[localIdx];
239 CornerStorage corners_;
240 GlobalPosition center_;
242 GridIndexType elementIndex_;
243 LocalIndexType vIdxLocal_;
244 LocalIndexType elemLocalScvIdx_;
245 GridIndexType dofIndex_;
248 LocalIndexType facetIdx_;
Define some often used mathematical functions.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for a sub control volume.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
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
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:50
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:54
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:53
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:51
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, ScvMLGTraits< Scalar > > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:80
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:77
typename ScvMLGTraits< Scalar >::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:81
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:82
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:78
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:79
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:58
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(dim)) > Type
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:65
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:71
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:73
static const bool v
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:72
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:96
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:182
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:198
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:194
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:219
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:214
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:231
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:186
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:206
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:111
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:210
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:148
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:190
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:118
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:202