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_; }
187 [[deprecated(
"Will be removed after release 3.9. Use periodicDofMap() instead.")]]
189 {
return periodicDofMap(); }
193 {
return { gg.cache_ }; }
197 class BoxGridGeometryCache
208 const BoxFVGridGeometry& gridGeometry()
const
209 {
return *gridGeometry_; }
212 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx)
const
213 {
return scvs_[eIdx]; }
216 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx)
const
217 {
return scvfs_[eIdx]; }
220 bool hasBoundaryScvf(GridIndexType eIdx)
const
221 {
return hasBoundaryScvf_[eIdx]; }
224 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx)
const
225 {
return scvfBoundaryGeometryKeys_.at(eIdx); }
232 hasBoundaryScvf_.clear();
233 scvfBoundaryGeometryKeys_.clear();
236 std::vector<std::vector<SubControlVolume>> scvs_;
237 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
238 std::vector<bool> hasBoundaryScvf_;
239 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
241 const BoxFVGridGeometry* gridGeometry_;
250 using GeometryHelper =
typename Cache::GeometryHelper;
256 const auto numElements = this->gridView().size(0);
257 cache_.scvs_.resize(numElements);
258 cache_.scvfs_.resize(numElements);
259 cache_.hasBoundaryScvf_.resize(numElements,
false);
261 boundaryDofIndices_.assign(numDofs(),
false);
265 numBoundaryScvf_ = 0;
267 for (
const auto& element : elements(this->gridView()))
270 const auto eIdx = this->elementMapper().index(element);
273 numScv_ += element.subEntities(dim);
274 numScvf_ += element.subEntities(dim-1);
277 auto elementGeometry = element.geometry();
278 const auto refElement = referenceElement(elementGeometry);
281 GeometryHelper geometryHelper(elementGeometry);
284 cache_.scvs_[eIdx].resize(elementGeometry.corners());
285 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
287 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
290 geometryHelper.getScvCorners(scvLocalIdx),
298 LocalIndexType scvfLocalIdx = 0;
299 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
300 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
303 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
304 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
306 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
309 geometryHelper.normal(corners, localScvIndices),
313 std::move(localScvIndices),
319 for (
const auto& intersection : intersections(this->gridView(), element))
321 if (intersection.boundary() && !intersection.neighbor())
323 const auto isGeometry = intersection.geometry();
324 cache_.hasBoundaryScvf_[eIdx] =
true;
327 numScvf_ += isGeometry.corners();
328 numBoundaryScvf_ += isGeometry.corners();
330 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
333 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
334 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
336 cache_.scvfs_[eIdx].emplace_back(
337 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
338 intersection.centerUnitOuterNormal(),
343 std::move(localScvIndices),
347 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
348 static_cast<LocalIndexType
>(intersection.indexInInside()),
349 static_cast<LocalIndexType
>(isScvfLocalIdx)
358 const auto fIdx = intersection.indexInInside();
359 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
360 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
362 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
363 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
364 boundaryDofIndices_[vIdxGlobal] =
true;
369 else if (intersection.boundary() && intersection.neighbor())
374 const auto fIdx = intersection.indexInInside();
375 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
376 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
377 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
379 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
380 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
381 const auto vPos = elementGeometry.corner(vIdx);
383 const auto& outside = intersection.outside();
384 const auto outsideGeometry = outside.geometry();
385 for (
const auto& isOutside : intersections(this->gridView(), outside))
388 if (isOutside.boundary() && isOutside.neighbor())
390 const auto fIdxOutside = isOutside.indexInInside();
391 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
392 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
394 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
395 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
396 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
397 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
398 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
408 if (this->isPeriodic() && this->gridView().comm().size() > 1)
409 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
412 const FeCache feCache_;
415 std::size_t numScvf_;
416 std::size_t numBoundaryScvf_;
419 std::vector<bool> boundaryDofIndices_;
422 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
434template<
class Scalar,
class GV,
class Traits>
442 static const int dim = GV::dimension;
443 static const int dimWorld = GV::dimensionworld;
445 using Element =
typename GV::template Codim<0>::Entity;
446 using CoordScalar =
typename GV::ctype;
466 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
486 {
return this->vertexMapper(); }
499 {
return numBoundaryScvf_; }
503 {
return this->vertexMapper().size(); }
509 ParentType::update(gridView);
516 ParentType::update(std::move(gridView));
526 {
return boundaryDofIndices_[dofIdx]; }
530 {
return periodicVertexMap_.count(dofIdx); }
534 {
return periodicVertexMap_.at(dofIdx); }
538 {
return periodicVertexMap_; }
541 [[deprecated(
"Will be removed after release 3.9. Use periodicDofMap() instead.")]]
543 {
return periodicDofMap(); }
547 {
return { gg.cache_ }; }
551 class BoxGridGeometryCache
562 const BoxFVGridGeometry& gridGeometry()
const
563 {
return *gridGeometry_; }
566 const BoxFVGridGeometry* gridGeometry_;
578 boundaryDofIndices_.assign(numDofs(),
false);
584 numBoundaryScvf_ = 0;
585 for (
const auto& element : elements(this->gridView()))
587 numScv_ += element.subEntities(dim);
588 numScvf_ += element.subEntities(dim-1);
590 const auto elementGeometry = element.geometry();
591 const auto refElement = referenceElement(elementGeometry);
594 for (
const auto& intersection : intersections(this->gridView(), element))
596 if (intersection.boundary() && !intersection.neighbor())
598 const auto isGeometry = intersection.geometry();
599 numScvf_ += isGeometry.corners();
600 numBoundaryScvf_ += isGeometry.corners();
604 const auto fIdx = intersection.indexInInside();
605 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
606 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
608 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
609 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
610 boundaryDofIndices_[vIdxGlobal] =
true;
615 else if (intersection.boundary() && intersection.neighbor())
620 const auto fIdx = intersection.indexInInside();
621 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
622 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
623 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
625 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
626 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
627 const auto vPos = elementGeometry.corner(vIdx);
629 const auto& outside = intersection.outside();
630 const auto outsideGeometry = outside.geometry();
631 for (
const auto& isOutside : intersections(this->gridView(), outside))
634 if (isOutside.boundary() && isOutside.neighbor())
636 const auto fIdxOutside = isOutside.indexInInside();
637 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
638 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
640 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
641 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
642 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
643 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
644 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
654 if (this->isPeriodic() && this->gridView().comm().size() > 1)
655 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
658 const FeCache feCache_;
663 std::size_t numScvf_;
664 std::size_t numBoundaryScvf_;
667 std::vector<bool> boundaryDofIndices_;
670 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:437
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:498
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:546
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:529
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:456
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:460
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:479
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:458
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:533
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:521
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:572
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:493
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:464
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:454
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:525
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:485
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:514
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:462
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:502
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:542
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:471
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:489
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:507
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:466
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:537
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:192
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:247
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:188
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
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:183
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