15#ifndef DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
16#define DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
21#include <dune/grid/common/mcmgmapper.hh>
22#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
38template<
class GV,
class T>
53template<
class Gr
idView>
60 template<
class Gr
idGeometry,
bool enableCache>
64 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
68 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
80 bool enableGridGeometryCache =
false,
91template<
class Scalar,
class GV,
class Traits>
100 using Element =
typename GV::template Codim<0>::Entity;
101 using CoordScalar =
typename GV::ctype;
102 static const int dim = GV::dimension;
103 static const int dimWorld = GV::dimensionworld;
121 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
126 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
128 const FacetGridView& facetGridView,
130 bool verbose =
false)
134 update_(facetGridView, codimOneGridAdapter, verbose);
139 {
return this->vertexMapper(); }
152 {
return numBoundaryScvf_; }
156 {
return this->vertexMapper().size(); }
171 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
173 const FacetGridView& facetGridView,
175 bool verbose =
false)
177 ParentType::update(gridView);
178 update_(facetGridView, codimOneGridAdapter, verbose);
182 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
184 const FacetGridView& facetGridView,
186 bool verbose =
false)
188 ParentType::update(std::move(gridView));
189 update_(facetGridView, codimOneGridAdapter, verbose);
198 {
return boundaryDofIndices_[dofIdx]; }
202 {
return interiorBoundaryDofIndices_[dofIdx]; }
210 { DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box facet coupling scheme"); }
214 {
return { gg.cache_ }; }
217 class BoxFacetCouplingGridGeometryCache
228 const BoxFacetCouplingFVGridGeometry& gridGeometry()
const
229 {
return *gridGeometry_; }
232 const std::vector<SubControlVolume>&
scvs(GridIndexType eIdx)
const
233 {
return scvs_[eIdx]; }
236 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx)
const
237 {
return scvfs_[eIdx]; }
247 std::vector<std::vector<SubControlVolume>> scvs_;
248 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
250 const BoxFacetCouplingFVGridGeometry* gridGeometry_;
256 using Cache = BoxFacetCouplingGridGeometryCache;
259 using GeometryHelper =
typename Cache::GeometryHelper;
261 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
262 void update_(
const FacetGridView& facetGridView,
264 bool verbose =
false)
267 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
272 const auto numDof = numDofs();
273 const auto numElements = this->gridView().size(0);
274 cache_.scvs_.resize(numElements);
275 cache_.scvfs_.resize(numElements);
276 boundaryDofIndices_.assign(numDof,
false);
277 interiorBoundaryDofIndices_.assign(numDof,
false);
282 numBoundaryScvf_ = 0;
283 for (
const auto& element : elements(this->gridView()))
285 auto eIdx = this->elementMapper().index(element);
288 numScv_ += element.subEntities(dim);
289 numScvf_ += element.subEntities(dim-1);
292 auto elementGeometry = element.geometry();
293 const auto refElement = referenceElement(elementGeometry);
296 GeometryHelper geometryHelper(elementGeometry);
299 cache_.scvs_[eIdx].resize(elementGeometry.corners());
300 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
301 cache_.scvs_[eIdx][scvLocalIdx] =
SubControlVolume(geometryHelper.getScvCorners(scvLocalIdx),
304 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
307 LocalIndexType scvfLocalIdx = 0;
308 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
309 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
312 std::array<LocalIndexType, 2> localScvIndices{{
313 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
314 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
321 std::move(localScvIndices));
326 std::vector<unsigned int> handledFacets;
327 for (
const auto& intersection : intersections(this->gridView(), element))
329 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
332 handledFacets.push_back(intersection.indexInInside());
335 const auto isGeometry = intersection.geometry();
336 const auto numFaceCorners = isGeometry.corners();
337 const auto idxInInside = intersection.indexInInside();
338 const auto boundary = intersection.boundary();
340 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
341 for (
int i = 0; i < numFaceCorners; ++i)
342 vIndicesLocal[i] =
static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, i, dim));
344 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
345 for (
int i = 0; i < numFaceCorners; ++i)
346 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
352 if (boundary && intersection.neighbor())
353 DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box facet coupling scheme");
355 if (isOnFacet || boundary)
358 numScvf_ += numFaceCorners;
359 numBoundaryScvf_ += int(boundary)*numFaceCorners;
361 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
364 std::array<LocalIndexType, 2> localScvIndices{{
365 vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]
369 cache_.scvfs_[eIdx].emplace_back(geometryHelper,
373 std::move(localScvIndices),
378 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
379 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
380 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
390 const FeCache feCache_;
394 std::size_t numScvf_;
395 std::size_t numBoundaryScvf_;
398 std::vector<bool> boundaryDofIndices_;
399 std::vector<bool> interiorBoundaryDofIndices_;
411template<
class Scalar,
class GV,
class Traits>
420 static const int dim = GV::dimension;
421 static const int dimWorld = GV::dimensionworld;
423 using Element =
typename GV::template Codim<0>::Entity;
424 using Intersection =
typename GV::Intersection;
425 using CoordScalar =
typename GV::ctype;
443 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
448 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
450 const FacetGridView& facetGridView,
452 bool verbose =
false)
454 , facetMapper_(gridView,
Dune::mcmgLayout(
Dune::template Codim<1>()))
457 update_(facetGridView, codimOneGridAdapter, verbose);
463 {
return this->vertexMapper(); }
476 {
return numBoundaryScvf_; }
480 {
return this->vertexMapper().size(); }
495 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
497 const FacetGridView& facetGridView,
499 bool verbose =
false)
501 ParentType::update(gridView);
502 updateFacetMapper_();
503 update_(facetGridView, codimOneGridAdapter, verbose);
507 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
509 const FacetGridView& facetGridView,
511 bool verbose =
false)
513 ParentType::update(std::move(gridView));
514 updateFacetMapper_();
515 update_(facetGridView, codimOneGridAdapter, verbose);
524 {
return boundaryDofIndices_[dofIdx]; }
528 {
return interiorBoundaryDofIndices_[dofIdx]; }
532 {
return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
540 { DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the facet coupling scheme"); }
544 {
return { gg.cache_ }; }
548 class BoxFacetCouplingGridGeometryCache
559 const BoxFacetCouplingFVGridGeometry& gridGeometry()
const
560 {
return *gridGeometry_; }
563 const BoxFacetCouplingFVGridGeometry* gridGeometry_;
569 using Cache = BoxFacetCouplingGridGeometryCache;
573 void updateFacetMapper_()
575 facetMapper_.update(this->gridView());
578 template<
class FacetGr
idView,
class CodimOneGr
idAdapter>
579 void update_(
const FacetGridView& facetGridView,
584 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
590 numBoundaryScvf_ = 0;
592 const auto numDof = numDofs();
593 boundaryDofIndices_.assign(numDof,
false);
594 interiorBoundaryDofIndices_.assign(numDof,
false);
595 facetIsOnInteriorBoundary_.assign(this->gridView().size(1),
false);
596 for (
const auto& element : elements(this->gridView()))
598 numScv_ += element.subEntities(dim);
599 numScvf_ += element.subEntities(dim-1);
601 const auto elementGeometry = element.geometry();
602 const auto refElement = referenceElement(elementGeometry);
606 std::vector<unsigned int> handledFacets;
607 for (
const auto& intersection : intersections(this->gridView(), element))
609 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
612 handledFacets.push_back(intersection.indexInInside());
615 const auto isGeometry = intersection.geometry();
616 const auto numFaceCorners = isGeometry.corners();
617 const auto idxInInside = intersection.indexInInside();
618 const auto boundary = intersection.boundary();
620 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
621 for (
int i = 0; i < numFaceCorners; ++i)
622 vIndicesLocal[i] =
static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, i, dim));
624 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
625 for (
int i = 0; i < numFaceCorners; ++i)
626 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
632 if (boundary && intersection.neighbor())
633 DUNE_THROW(Dune::InvalidStateException,
"Periodic boundaries are not supported by the box facet coupling scheme");
635 if (isOnFacet || boundary)
637 numScvf_ += numFaceCorners;
638 numBoundaryScvf_ += int(boundary)*numFaceCorners;
641 for (
int i = 0; i < numFaceCorners; ++i)
643 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
644 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
647 interiorBoundaryDofIndices_[ dofIndex ] =
true;
648 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] =
true;
656 const FeCache feCache_;
661 std::size_t numScvf_;
662 std::size_t numBoundaryScvf_;
665 std::vector<bool> boundaryDofIndices_;
666 std::vector<bool> interiorBoundaryDofIndices_;
669 typename Traits::FacetMapper facetMapper_;
670 std::vector<bool> facetIsOnInteriorBoundary_;
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
Base class for the element-local finite volume geometry for box models in the context of models consi...
Definition: multidomain/facet/box/fvelementgeometry.hh:36
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:414
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:435
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:449
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:466
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:470
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:479
void update(GridView &&gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:508
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:475
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:443
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:519
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:441
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:462
friend LocalView localView(const BoxFacetCouplingFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: multidomain/facet/box/fvgridgeometry.hh:543
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: multidomain/facet/box/fvgridgeometry.hh:433
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:539
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:535
BoxFacetCouplingGridGeometryCache Cache
Definition: multidomain/facet/box/fvgridgeometry.hh:569
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:437
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:439
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:523
void update(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:496
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:527
bool isOnInteriorBoundary(const Element &element, const Intersection &intersection) const
returns true if an intersection is on an interior boundary
Definition: multidomain/facet/box/fvgridgeometry.hh:531
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:94
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:201
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:209
void update(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:172
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:121
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:117
BoxFacetCouplingGridGeometryCache Cache
Definition: multidomain/facet/box/fvgridgeometry.hh:256
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:142
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:205
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:197
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:119
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:151
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:138
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:113
friend LocalView localView(const BoxFacetCouplingFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: multidomain/facet/box/fvgridgeometry.hh:213
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:146
void update(GridView &&gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:183
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:193
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:155
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:127
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: multidomain/facet/box/fvgridgeometry.hh:111
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:115
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:82
Class for a sub control volume face in the box method, i.e a part of the boundary of a sub control vo...
Definition: multidomain/facet/box/subcontrolvolumeface.hh:39
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:257
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:58
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:40
bool composeFacetElement(const IndexStorage &bulkVertexIndices) const
Returns true if a given set of bulk vertex indices make up a facet grid element.
Definition: codimonegridadapter.hh:180
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...
Definition: vertexmapper.hh:126
the sub control volume for the box scheme
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for a sub control volume face of the box method in the context of of models considering co...
Dune::Std::detected_or_t< Dumux::BoxGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxFacetCouplingGeometryHelper_t
Definition: multidomain/facet/box/fvgridgeometry.hh:43
typename T::GeometryHelper SpecifiesGeometryHelper
Definition: basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
Definition: common/pdesolver.hh:24
The default traits for the finite volume grid geometry of the box scheme with coupling occurring acro...
Definition: multidomain/facet/box/fvgridgeometry.hh:55
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:68
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:64
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26