14#ifndef DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
18#include <unordered_map>
22#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
37template<
class GV,
class T>
51template<
class Gr
idView,
class MapperTraits = DefaultMapperTraits<Gr
idView>>
58 template<
class Gr
idGeometry,
bool enableCache>
70 bool enableGridGeometryCache =
false,
80template<
class Scalar,
class GV,
class Traits>
89 using Element =
typename GV::template Codim<0>::Entity;
90 using CoordScalar =
typename GV::ctype;
91 static const int dim = GV::dimension;
92 static const int dimWorld = GV::dimensionworld;
112 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
132 {
return this->vertexMapper(); }
145 {
return numBoundaryScvf_; }
149 {
return this->vertexMapper().size(); }
155 ParentType::update(gridView);
162 ParentType::update(std::move(gridView));
172 {
return boundaryDofIndices_[dofIdx]; }
176 {
return periodicVertexMap_.count(dofIdx); }
180 {
return periodicVertexMap_.at(dofIdx); }
184 {
return periodicVertexMap_; }
188 {
return { gg.cache_ }; }
192 class BoxGridGeometryCache
203 const BoxFVGridGeometry& gridGeometry()
const
204 {
return *gridGeometry_; }
207 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx)
const
208 {
return scvs_[eIdx]; }
211 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx)
const
212 {
return scvfs_[eIdx]; }
215 bool hasBoundaryScvf(GridIndexType eIdx)
const
216 {
return hasBoundaryScvf_[eIdx]; }
219 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx)
const
220 {
return scvfBoundaryGeometryKeys_.at(eIdx); }
227 hasBoundaryScvf_.clear();
228 scvfBoundaryGeometryKeys_.clear();
231 std::vector<std::vector<SubControlVolume>> scvs_;
232 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
233 std::vector<bool> hasBoundaryScvf_;
234 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
236 const BoxFVGridGeometry* gridGeometry_;
245 using GeometryHelper =
typename Cache::GeometryHelper;
251 const auto numElements = this->gridView().size(0);
252 cache_.scvs_.resize(numElements);
253 cache_.scvfs_.resize(numElements);
254 cache_.hasBoundaryScvf_.resize(numElements,
false);
256 boundaryDofIndices_.assign(numDofs(),
false);
260 numBoundaryScvf_ = 0;
262 for (
const auto& element : elements(this->gridView()))
265 const auto eIdx = this->elementMapper().index(element);
268 numScv_ += element.subEntities(dim);
269 numScvf_ += element.subEntities(dim-1);
272 auto elementGeometry = element.geometry();
273 const auto refElement = referenceElement(elementGeometry);
276 GeometryHelper geometryHelper(elementGeometry);
279 cache_.scvs_[eIdx].resize(elementGeometry.corners());
280 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
282 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
285 geometryHelper.getScvCorners(scvLocalIdx),
293 LocalIndexType scvfLocalIdx = 0;
294 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
295 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
298 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
299 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
301 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
304 geometryHelper.normal(corners, localScvIndices),
308 std::move(localScvIndices),
314 for (
const auto& intersection : intersections(this->gridView(), element))
316 if (intersection.boundary() && !intersection.neighbor())
318 const auto isGeometry = intersection.geometry();
319 cache_.hasBoundaryScvf_[eIdx] =
true;
322 numScvf_ += isGeometry.corners();
323 numBoundaryScvf_ += isGeometry.corners();
325 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
328 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
329 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
331 cache_.scvfs_[eIdx].emplace_back(
332 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
333 intersection.centerUnitOuterNormal(),
338 std::move(localScvIndices),
342 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
343 static_cast<LocalIndexType
>(intersection.indexInInside()),
344 static_cast<LocalIndexType
>(isScvfLocalIdx)
353 const auto fIdx = intersection.indexInInside();
354 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
355 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
357 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
358 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
359 boundaryDofIndices_[vIdxGlobal] =
true;
364 else if (intersection.boundary() && intersection.neighbor())
369 const auto fIdx = intersection.indexInInside();
370 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
371 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
372 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
374 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
375 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
376 const auto vPos = elementGeometry.corner(vIdx);
378 const auto& outside = intersection.outside();
379 const auto outsideGeometry = outside.geometry();
380 for (
const auto& isOutside : intersections(this->gridView(), outside))
383 if (isOutside.boundary() && isOutside.neighbor())
385 const auto fIdxOutside = isOutside.indexInInside();
386 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
387 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
389 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
390 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
391 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
392 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
393 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
403 if (this->isPeriodic() && this->gridView().comm().size() > 1)
404 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
407 const FeCache feCache_;
410 std::size_t numScvf_;
411 std::size_t numBoundaryScvf_;
414 std::vector<bool> boundaryDofIndices_;
417 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
429template<
class Scalar,
class GV,
class Traits>
437 static const int dim = GV::dimension;
438 static const int dimWorld = GV::dimensionworld;
440 using Element =
typename GV::template Codim<0>::Entity;
441 using CoordScalar =
typename GV::ctype;
461 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
481 {
return this->vertexMapper(); }
494 {
return numBoundaryScvf_; }
498 {
return this->vertexMapper().size(); }
504 ParentType::update(gridView);
511 ParentType::update(std::move(gridView));
521 {
return boundaryDofIndices_[dofIdx]; }
525 {
return periodicVertexMap_.count(dofIdx); }
529 {
return periodicVertexMap_.at(dofIdx); }
533 {
return periodicVertexMap_; }
537 {
return { gg.cache_ }; }
541 class BoxGridGeometryCache
552 const BoxFVGridGeometry& gridGeometry()
const
553 {
return *gridGeometry_; }
556 const BoxFVGridGeometry* gridGeometry_;
568 boundaryDofIndices_.assign(numDofs(),
false);
574 numBoundaryScvf_ = 0;
575 for (
const auto& element : elements(this->gridView()))
577 numScv_ += element.subEntities(dim);
578 numScvf_ += element.subEntities(dim-1);
580 const auto elementGeometry = element.geometry();
581 const auto refElement = referenceElement(elementGeometry);
584 for (
const auto& intersection : intersections(this->gridView(), element))
586 if (intersection.boundary() && !intersection.neighbor())
588 const auto isGeometry = intersection.geometry();
589 numScvf_ += isGeometry.corners();
590 numBoundaryScvf_ += isGeometry.corners();
594 const auto fIdx = intersection.indexInInside();
595 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
596 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
598 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
599 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
600 boundaryDofIndices_[vIdxGlobal] =
true;
605 else if (intersection.boundary() && intersection.neighbor())
610 const auto fIdx = intersection.indexInInside();
611 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
612 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
613 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
615 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
616 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
617 const auto vPos = elementGeometry.corner(vIdx);
619 const auto& outside = intersection.outside();
620 const auto outsideGeometry = outside.geometry();
621 for (
const auto& isOutside : intersections(this->gridView(), outside))
624 if (isOutside.boundary() && isOutside.neighbor())
626 const auto fIdxOutside = isOutside.indexInInside();
627 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
628 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
630 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
631 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
632 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
633 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
634 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
644 if (this->isPeriodic() && this->gridView().comm().size() > 1)
645 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
648 const FeCache feCache_;
653 std::size_t numScvf_;
654 std::size_t numBoundaryScvf_;
657 std::vector<bool> boundaryDofIndices_;
660 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
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 finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:41
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:432
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:493
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:536
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:524
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:451
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:455
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:474
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:453
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:528
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:516
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:562
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:488
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:459
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:449
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:520
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:480
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:509
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:457
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:497
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:532
BoxFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg)
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/box/fvgridgeometry.hh:466
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:484
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:502
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:461
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:83
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:187
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:242
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:183
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:153
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:179
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:148
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:135
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:100
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:167
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:144
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:112
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:106
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:131
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:125
BoxFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg)
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/box/fvgridgeometry.hh:117
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:108
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:102
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:171
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:139
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:175
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:160
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:104
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:110
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:72
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:257
Class for a sub control volume face in the box method, i.e a part of the boundary of a sub control vo...
Definition: discretization/box/subcontrolvolumeface.hh:62
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:60
Defines the default element and vertex mapper types.
Base class for the local finite volume geometry for box models This builds up the sub control volumes...
the sub control volume for the box scheme
Base class for a sub control volume face.
Helper classes to compute the integration elements.
Dune::Std::detected_or_t< Dumux::BasicGridGeometry< GV, typename T::ElementMapper, typename T::VertexMapper >, Detail::SpecifiesBaseGridGeometry, T > BasicGridGeometry_t
Type of the basic grid geometry implementation used as backend.
Definition: basegridgeometry.hh:42
The available discretization methods in Dumux.
Dune::Std::detected_or_t< Dumux::BoxGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxGeometryHelper_t
Definition: discretization/box/fvgridgeometry.hh:42
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
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition: discretization/box/fvgridgeometry.hh:54
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26