25#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
28#include <dune/common/reservedvector.hh>
29#include <dune/geometry/multilineargeometry.hh>
47template<
class Gr
idView>
50 using Grid =
typename GridView::Grid;
52 static const int dim = Grid::dimension;
53 static const int dimWorld = Grid::dimensionworld;
61 template<
int mydim,
int cdim >
64 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(
dim)) >;
71 static const bool v =
true;
72 static const unsigned int topologyId = Dune::Impl::CubeTopology< mydim >::type::id;
79 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, ScvMLGTraits<Scalar>>;
98 using Geometry =
typename T::Geometry;
99 using GridIndexType =
typename T::GridIndexType;
100 using LocalIndexType =
typename T::LocalIndexType;
101 using Scalar =
typename T::Scalar;
102 using GlobalPosition =
typename T::GlobalPosition;
103 using CornerStorage =
typename T::CornerStorage;
104 enum { dim = Geometry::mydimension };
106 static_assert(dim == 2 || dim == 3,
"Box-Dfm sub-control volume only implemented in 2d or 3d");
116 template<
class GeometryHelper>
118 LocalIndexType scvIdx,
121 : isFractureScv_(false)
122 , corners_(geometryHelper.getScvCorners(scvIdx))
124 , volume_(geometryHelper.scvVolume(corners_))
127 , elemLocalScvIdx_(scvIdx)
132 for (
const auto&
corner : corners_)
134 center_ /= corners_.size();
146 template<
class GeometryHelper,
class Intersection>
148 const Intersection& intersection,
149 const typename Intersection::Geometry& isGeometry,
150 LocalIndexType indexInIntersection,
151 LocalIndexType vIdxLocal,
152 LocalIndexType elemLocalScvIdx,
153 LocalIndexType elemLocalFacetIdx,
156 : isFractureScv_(true)
161 , vIdxLocal_(vIdxLocal)
162 , elemLocalScvIdx_(elemLocalScvIdx)
164 , facetIdx_(elemLocalFacetIdx)
167 auto corners = geometryHelper.getBoundaryScvfCorners(intersection, isGeometry, indexInIntersection);
168 const auto numCorners = corners.size();
169 corners_.resize(numCorners);
170 for (
unsigned int i = 0; i < numCorners; ++i)
171 corners_[i] = corners[i];
174 volume_ = geometryHelper.scvfArea(corners);
175 for (
const auto&
corner : corners_)
177 center_ /= corners_.size();
190 {
return vIdxLocal_; }
194 {
return elemLocalScvIdx_; }
198 { assert(isFractureScv_);
return facetIdx_; }
202 {
return dofIndex_; }
206 {
return corners_[0]; }
210 {
return elementIndex_; }
214 {
return isFractureScv_; }
221 DUNE_THROW(Dune::InvalidStateException,
"The geometry object cannot be defined for fracture scvs "
222 "because the number of known corners is insufficient. "
223 "You can do this manually by extract the corners from this scv "
224 "and extrude them by the corresponding aperture. ");
226 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
230 const GlobalPosition&
corner(LocalIndexType localIdx)
const
232 assert(localIdx < corners_.size() &&
"provided index exceeds the number of corners");
233 return corners_[localIdx];
238 CornerStorage corners_;
239 GlobalPosition center_;
241 GridIndexType elementIndex_;
242 LocalIndexType vIdxLocal_;
243 LocalIndexType elemLocalScvIdx_;
244 GridIndexType dofIndex_;
247 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.
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:49
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:53
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:52
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:50
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, ScvMLGTraits< Scalar > > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:79
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:76
typename ScvMLGTraits< Scalar >::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:80
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:81
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:77
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:78
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:57
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:63
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(dim)) > Type
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:70
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:72
static const bool v
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:71
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:95
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:181
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:197
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:193
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:218
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:213
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:230
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:185
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:205
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:110
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:209
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:147
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:189
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:117
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:201