29#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
30#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
32#include <unordered_map>
34#include <dune/geometry/referenceelements.hh>
35#include <dune/localfunctions/lagrange/pqkfactory.hh>
36#include <dune/geometry/multilineargeometry.hh>
37#include <dune/grid/common/mcmgmapper.hh>
59template<
class Gr
idView,
class MapperTraits = DefaultMapperTraits<Gr
idView>>
66 template<
class Gr
idGeometry,
bool enableCache>
70 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
83 bool enableGridGeometryCache =
false,
97template<
class Scalar,
class GV,
class Traits>
103 using GridIndexType =
typename GV::IndexSet::IndexType;
105 using Element =
typename GV::template Codim<0>::Entity;
106 using CoordScalar =
typename GV::ctype;
107 static const int dim = GV::dimension;
108 static const int dimWorld = GV::dimensionworld;
109 static_assert(dim == 2 || dim == 3,
"The box-dfm GridGeometry is only implemented in 2 or 3 dimensions.");
111 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
112 using FaceReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim-1>;
115 typename Traits::SubControlVolume,
116 typename Traits::SubControlVolumeFace>;
131 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
155 {
return numBoundaryScvf_; }
159 {
return this->
gridView().size(dim); }
162 template<
class FractureGr
idAdapter >
163 void update(
const FractureGridAdapter& fractureGridAdapter)
170 auto numElements = this->
gridView().size(0);
171 scvs_.resize(numElements);
172 scvfs_.resize(numElements);
174 boundaryDofIndices_.assign(
numDofs(),
false);
175 fractureDofIndices_.assign(this->
gridView.size(dim),
false);
179 numBoundaryScvf_ = 0;
187 numScv_ +=
element.subEntities(dim);
188 numScvf_ +=
element.subEntities(dim-1);
191 auto elementGeometry =
element.geometry();
192 const auto referenceElement = ReferenceElements::general(elementGeometry.type());
195 GeometryHelper geometryHelper(elementGeometry);
198 scvs_[eIdx].resize(elementGeometry.corners());
199 using LocalIndexType =
typename SubControlVolumeFace::Traits::LocalIndexType;
200 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
211 LocalIndexType scvfLocalIdx = 0;
212 scvfs_[eIdx].resize(
element.subEntities(dim-1));
213 for (; scvfLocalIdx <
element.subEntities(dim-1); ++scvfLocalIdx)
216 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(referenceElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
217 static_cast<LocalIndexType
>(referenceElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
223 std::move(localScvIndices));
237 LocalIndexType scvLocalIdx =
element.subEntities(dim);
238 for (
const auto& intersection : intersections(this->
gridView(),
element))
241 const auto& isGeometry = intersection.geometry();
242 const auto numCorners = isGeometry.corners();
243 const auto idxInInside = intersection.indexInInside();
245 std::vector<GridIndexType> isVertexIndices(numCorners);
246 for (
unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
248 referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim),
251 if (intersection.boundary() && !intersection.neighbor())
253 numScvf_ += isGeometry.corners();
254 numBoundaryScvf_ += isGeometry.corners();
256 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
259 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(referenceElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
260 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
261 scvfs_[eIdx].emplace_back(geometryHelper,
266 std::move(localScvIndices));
271 const auto numFaceVerts = referenceElement.size(idxInInside, 1, dim);
272 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
274 const auto vIdx = referenceElement.subEntity(idxInInside, 1, localVIdx, dim);
276 boundaryDofIndices_[vIdxGlobal] =
true;
280 else if (intersection.boundary() && intersection.neighbor())
281 DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box-dfm scheme");
284 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
286 for (
auto vIdx : isVertexIndices)
287 fractureDofIndices_[vIdx] =
true;
290 numScv_ += numCorners;
291 const auto curNumScvs = scvs_[eIdx].size();
292 scvs_[eIdx].reserve(curNumScvs+numCorners);
293 for (
unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
294 scvs_[eIdx].emplace_back(geometryHelper,
298 static_cast<LocalIndexType
>(referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
302 isVertexIndices[vIdxLocal]);
307 const auto& faceReferenceElement = FaceReferenceElements::general(isGeometry.type());
308 for (
unsigned int edgeIdx = 0; edgeIdx < faceReferenceElement.size(1); ++edgeIdx)
311 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(faceReferenceElement.subEntity(edgeIdx, 1, 0, dim-1)),
312 static_cast<LocalIndexType
>(faceReferenceElement.subEntity(edgeIdx, 1, 1, dim-1))});
315 std::for_each( localScvIndices.begin(),
316 localScvIndices.end(),
317 [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
321 scvfs_[eIdx].emplace_back(geometryHelper,
326 std::move(localScvIndices),
327 intersection.boundary());
335 std::vector<LocalIndexType> localScvIndices({0, 1});
338 std::for_each( localScvIndices.begin(),
339 localScvIndices.end(),
340 [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
344 scvfs_[eIdx].emplace_back(geometryHelper,
349 std::move(localScvIndices),
350 intersection.boundary());
360 const std::vector<SubControlVolume>&
scvs(GridIndexType eIdx)
const {
return scvs_[eIdx]; }
362 const std::vector<SubControlVolumeFace>&
scvfs(GridIndexType eIdx)
const {
return scvfs_[eIdx]; }
364 bool dofOnBoundary(
unsigned int dofIdx)
const {
return boundaryDofIndices_[dofIdx]; }
366 bool dofOnFracture(
unsigned int dofIdx)
const {
return fractureDofIndices_[dofIdx]; }
372 { DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box-dfm scheme"); }
376 {
return std::unordered_map<std::size_t, std::size_t>(); }
379 const FeCache feCache_;
381 std::vector<std::vector<SubControlVolume>> scvs_;
382 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
386 std::size_t numScvf_;
387 std::size_t numBoundaryScvf_;
390 std::vector<bool> boundaryDofIndices_;
391 std::vector<bool> fractureDofIndices_;
401template<
class Scalar,
class GV,
class Traits>
407 using GridIndexType =
typename GV::IndexSet::IndexType;
409 static const int dim = GV::dimension;
410 static const int dimWorld = GV::dimensionworld;
412 using Element =
typename GV::template Codim<0>::Entity;
413 using Intersection =
typename GV::Intersection;
414 using CoordScalar =
typename GV::ctype;
416 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
417 using FaceReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim-1>;
432 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
458 {
return numBoundaryScvf_; }
462 {
return this->
gridView().size(dim); }
465 template<
class FractureGr
idAdapter >
466 void update(
const FractureGridAdapter& fractureGridAdapter)
470 boundaryDofIndices_.assign(
numDofs(),
false);
471 fractureDofIndices_.assign(
numDofs(),
false);
472 facetOnFracture_.assign(this->
gridView().size(1),
false);
478 numBoundaryScvf_ = 0;
481 numScv_ +=
element.subEntities(dim);
482 numScvf_ +=
element.subEntities(dim-1);
484 const auto elementGeometry =
element.geometry();
485 const auto referenceElement = ReferenceElements::general(elementGeometry.type());
488 for (
const auto& intersection : intersections(this->
gridView(),
element))
491 const auto& isGeometry = intersection.geometry();
492 const auto numCorners = isGeometry.corners();
493 const auto idxInInside = intersection.indexInInside();
495 std::vector<GridIndexType> isVertexIndices(numCorners);
496 for (
unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
498 referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim),
501 if (intersection.boundary() && !intersection.neighbor())
503 numScvf_ += numCorners;
504 numBoundaryScvf_ += numCorners;
508 const auto fIdx = intersection.indexInInside();
509 const auto numFaceVerts = referenceElement.size(fIdx, 1, dim);
510 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
512 const auto vIdx = referenceElement.subEntity(fIdx, 1, localVIdx, dim);
514 boundaryDofIndices_[vIdxGlobal] =
true;
519 else if (intersection.boundary() && intersection.neighbor())
520 DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box-dfm scheme");
523 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
525 facetOnFracture_[facetMapper_.subIndex(
element, idxInInside, 1)] =
true;
526 for (
auto vIdx : isVertexIndices)
527 fractureDofIndices_[vIdx] =
true;
529 const auto isGeometry = intersection.geometry();
530 numScv_ += isGeometry.corners();
531 numScvf_ += dim == 3 ? FaceReferenceElements::general(isGeometry.type()).size(1) : 1;
540 bool dofOnBoundary(
unsigned int dofIdx)
const {
return boundaryDofIndices_[dofIdx]; }
542 bool dofOnFracture(
unsigned int dofIdx)
const {
return fractureDofIndices_[dofIdx]; }
548 {
return facetOnFracture_[facetMapper_.subIndex(
element, intersection.indexInInside(), 1)]; }
552 { DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box-dfm scheme"); }
556 {
return std::unordered_map<std::size_t, std::size_t>(); }
560 const FeCache feCache_;
565 std::size_t numScvf_;
566 std::size_t numBoundaryScvf_;
569 std::vector<bool> boundaryDofIndices_;
570 std::vector<bool> fractureDofIndices_;
573 typename Traits::FacetMapper facetMapper_;
574 std::vector<bool> facetOnFracture_;
Defines the default element and vertex mapper types.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
The available discretization methods in Dumux.
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
BaseGridGeometry(const GridView &gridView)
Constructor computes the bouding box of the entire domain, for e.g. setting boundary conditions.
Definition basegridgeometry.hh:77
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ box
Definition method.hh:38
Definition common/pdesolver.hh:35
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices for constant grids.
Definition basegridgeometry.hh:119
Element element(GridIndexType eIdx) const
Get an element from a global element index.
Definition basegridgeometry.hh:163
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices for constant grids.
Definition basegridgeometry.hh:113
const GridView & gridView() const
Return the gridView this grid geometry object lives on.
Definition basegridgeometry.hh:107
void update()
Update all fvElementGeometries (do this again after grid adaption).
Definition basegridgeometry.hh:90
Base class for the finite volume geometry vector for box discrete fracture model.
Definition porousmediumflow/boxdfm/fvelementgeometry.hh:52
The default traits for the box finite volume grid geometry.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:62
BoxDfmSubControlVolume< GridView > SubControlVolume
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:63
BoxDfmSubControlVolumeFace< GridView > SubControlVolumeFace
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:64
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:70
BoxDfmFVElementGeometry< GridGeometry, enableCache > LocalView
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:67
Base class for the finite volume geometry vector for box schemes.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:85
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:125
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:375
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:362
std::size_t numDofs() const
The total number of degrees of freedom.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:158
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:129
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:136
GV GridView
Export the grid view type.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:133
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:131
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:360
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type).
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:123
const DofMapper & dofMapper() const
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:141
std::size_t numScvf() const
The total number of sun control volume faces.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:149
std::size_t numScv() const
The total number of sub control volumes.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:145
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:371
std::size_t numBoundaryScvf() const
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:154
static constexpr DiscretizationMethod discMethod
Export discretization method.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:120
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:364
void update(const FractureGridAdapter &fractureGridAdapter)
Update all fvElementGeometries (do this again after grid adaption).
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:163
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:358
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:127
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:368
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:366
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:540
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:544
std::size_t numScv() const
The total number of sub control volumes.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:448
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:547
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:555
std::size_t numDofs() const
The total number of degrees of freedom.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:461
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:437
static constexpr DiscretizationMethod discMethod
export discretization method
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:421
std::size_t numBoundaryScvf() const
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:457
GV GridView
export the grid view type
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:434
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:430
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:424
std::size_t numScvf() const
The total number of sun control volume faces.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:452
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:432
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:538
const DofMapper & dofMapper() const
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:444
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:428
void update(const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (do this again after grid adaption)
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:466
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:426
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:542
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:551
Create sub control volumes and sub control volume face geometries.
Definition geometryhelper.hh:35
the sub control volume for the box discrete fracture scheme
Definition porousmediumflow/boxdfm/subcontrolvolume.hh:96
Class for a sub control volume face in the box discrete fracture method, i.e a part of the boundary o...
Definition porousmediumflow/boxdfm/subcontrolvolumeface.hh:100
Base class for the local finite volume geometry for the box discrete fracture model.
the sub control volume for the box discrete fracture scheme
The sub control volume face class for the box discrete fracture model.