14#ifndef DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
18#include <unordered_map>
23#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
41template<
class GV,
class T>
53template<
class Gr
idView,
class ScvRule = Dumux::QuadratureRules::M
idpo
intQuadrature,
class ScvfRule = Dumux::QuadratureRules::M
idpo
intQuadrature>
62template<
class Gr
idView,
class MapperTraits = DefaultMapperTraits<Gr
idView>,
class QuadratureTraits = BoxQuadratureTraits<Gr
idView>>
64:
public MapperTraits,
public QuadratureTraits
69 template<
class Gr
idGeometry,
bool enableCache>
81 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;
125 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
137 : ParentType(std::move(gg))
139 , periodicGridTraits_(this->
gridView().grid())
165 {
return numBoundaryScvf_; }
192 {
return boundaryDofIndices_[dofIdx]; }
196 {
return periodicDofMap_.count(dofIdx); }
200 {
return periodicDofMap_.at(dofIdx); }
204 {
return periodicDofMap_; }
208 {
return { gg.cache_ }; }
212 class BoxGridGeometryCache
223 const BoxFVGridGeometry& gridGeometry()
const
224 {
return *gridGeometry_; }
227 const std::vector<SubControlVolume>&
scvs(GridIndexType eIdx)
const
228 {
return scvs_[eIdx]; }
231 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx)
const
232 {
return scvfs_[eIdx]; }
235 bool hasBoundaryScvf(GridIndexType eIdx)
const
236 {
return hasBoundaryScvf_[eIdx]; }
239 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx)
const
240 {
return scvfBoundaryGeometryKeys_.at(eIdx); }
243 auto boundaryFaces(GridIndexType eIdx)
const -> std::span<const BoundaryFace>
245 if (
auto it = boundaryFaces_.find(eIdx); it != boundaryFaces_.end())
251 const auto& boundaryFaceScvfRanges(GridIndexType eIdx)
const
252 {
return boundaryFaceScvfRanges_.at(eIdx); }
259 hasBoundaryScvf_.clear();
260 scvfBoundaryGeometryKeys_.clear();
261 boundaryFaces_.clear();
262 boundaryFaceScvfRanges_.clear();
265 std::vector<std::vector<SubControlVolume>> scvs_;
266 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
267 std::vector<bool> hasBoundaryScvf_;
268 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
269 std::unordered_map<GridIndexType, Dune::ReservedVector<typename BoxFVGridGeometry::BoundaryFace, 2*dim>> boundaryFaces_;
270 std::unordered_map<GridIndexType, Dune::ReservedVector<std::array<LocalIndexType, 2>, 2*dim>> boundaryFaceScvfRanges_;
272 const BoxFVGridGeometry* gridGeometry_;
281 using GeometryHelper =
typename Cache::GeometryHelper;
287 const auto numElements = this->
gridView().size(0);
288 cache_.scvs_.resize(numElements);
289 cache_.scvfs_.resize(numElements);
290 cache_.hasBoundaryScvf_.resize(numElements,
false);
292 boundaryDofIndices_.assign(
numDofs(),
false);
296 numBoundaryScvf_ = 0;
304 numScv_ +=
element.subEntities(dim);
305 numScvf_ +=
element.subEntities(dim-1);
308 auto elementGeometry =
element.geometry();
309 const auto refElement = referenceElement(elementGeometry);
312 GeometryHelper geometryHelper(elementGeometry);
315 cache_.scvs_[eIdx].resize(elementGeometry.corners());
316 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
321 geometryHelper.getScvCorners(scvLocalIdx),
329 LocalIndexType scvfLocalIdx = 0;
330 cache_.scvfs_[eIdx].resize(
element.subEntities(dim-1));
331 for (; scvfLocalIdx <
element.subEntities(dim-1); ++scvfLocalIdx)
334 std::array<LocalIndexType, 2> localScvIndices{{
335 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
336 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
339 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
342 geometryHelper.normal(corners, localScvIndices),
345 std::move(localScvIndices)
350 LocalIndexType numBoundaryFaces = 0;
351 for (
const auto& intersection : intersections(this->
gridView(),
element))
353 if (intersection.boundary() && !intersection.neighbor())
355 const auto isGeometry = intersection.geometry();
356 cache_.hasBoundaryScvf_[eIdx] =
true;
362 intersection.centerUnitOuterNormal(),
364 static_cast<LocalIndexType
>(intersection.indexInInside()),
365 typename BoundaryFace::Traits::BoundaryFlag{intersection}
369 numScvf_ += isGeometry.corners();
370 numBoundaryScvf_ += isGeometry.corners();
373 cache_.boundaryFaceScvfRanges_[eIdx].push_back(std::array<LocalIndexType, 2>{{
375 static_cast<LocalIndexType
>(isGeometry.corners())
378 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
381 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
382 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
384 cache_.scvfs_[eIdx].emplace_back(
385 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
386 intersection.centerUnitOuterNormal(),
390 std::move(localScvIndices)
393 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
394 static_cast<LocalIndexType
>(intersection.indexInInside()),
395 static_cast<LocalIndexType
>(isScvfLocalIdx)
404 const auto fIdx = intersection.indexInInside();
405 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
406 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
408 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
410 boundaryDofIndices_[vIdxGlobal] =
true;
415 else if (periodicGridTraits_.isPeriodic(intersection))
420 const auto fIdx = intersection.indexInInside();
421 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
422 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
423 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
425 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
427 const auto vPos = elementGeometry.corner(vIdx);
429 const auto& outside = intersection.outside();
430 const auto outsideGeometry = outside.geometry();
431 for (
const auto& isOutside : intersections(this->
gridView(), outside))
434 if (periodicGridTraits_.isPeriodic(isOutside))
436 const auto fIdxOutside = isOutside.indexInInside();
437 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
438 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
440 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
441 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
442 const auto shift = std::abs((this->
bBoxMax()-this->
bBoxMin())*intersection.centerUnitOuterNormal());
443 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
444 periodicDofMap_[vIdxGlobal] = this->
vertexMapper().subIndex(outside, vIdxOutside, dim);
455 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
458 const FeCache feCache_;
461 std::size_t numScvf_;
462 std::size_t numBoundaryScvf_;
465 std::vector<bool> boundaryDofIndices_;
468 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
482template<
class Scalar,
class GV,
class Traits>
490 static const int dim = GV::dimension;
491 static const int dimWorld = GV::dimensionworld;
493 using Element =
typename GV::template Codim<0>::Entity;
494 using CoordScalar =
typename GV::ctype;
516 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
528 : ParentType(std::move(gg))
530 , periodicGridTraits_(this->
gridView().grid())
556 {
return numBoundaryScvf_; }
583 {
return boundaryDofIndices_[dofIdx]; }
587 {
return periodicDofMap_.count(dofIdx); }
591 {
return periodicDofMap_.at(dofIdx); }
595 {
return periodicDofMap_; }
599 {
return { gg.cache_ }; }
603 class BoxGridGeometryCache
614 const BoxFVGridGeometry& gridGeometry()
const
615 {
return *gridGeometry_; }
618 const BoxFVGridGeometry* gridGeometry_;
630 boundaryDofIndices_.assign(
numDofs(),
false);
636 numBoundaryScvf_ = 0;
639 numScv_ +=
element.subEntities(dim);
640 numScvf_ +=
element.subEntities(dim-1);
642 const auto elementGeometry =
element.geometry();
643 const auto refElement = referenceElement(elementGeometry);
646 for (
const auto& intersection : intersections(this->
gridView(),
element))
648 if (intersection.boundary() && !intersection.neighbor())
650 const auto isGeometry = intersection.geometry();
651 numScvf_ += isGeometry.corners();
652 numBoundaryScvf_ += isGeometry.corners();
656 const auto fIdx = intersection.indexInInside();
657 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
658 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
660 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
662 boundaryDofIndices_[vIdxGlobal] =
true;
667 else if (periodicGridTraits_.isPeriodic(intersection))
672 const auto fIdx = intersection.indexInInside();
673 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
674 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
675 for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
677 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
679 const auto vPos = elementGeometry.corner(vIdx);
681 const auto& outside = intersection.outside();
682 const auto outsideGeometry = outside.geometry();
683 for (
const auto& isOutside : intersections(this->
gridView(), outside))
686 if (periodicGridTraits_.isPeriodic(isOutside))
688 const auto fIdxOutside = isOutside.indexInInside();
689 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
690 for (
int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
692 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
693 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
694 const auto shift = std::abs((this->
bBoxMax()-this->
bBoxMin())*intersection.centerUnitOuterNormal());
695 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
696 periodicDofMap_[vIdxGlobal] = this->
vertexMapper().subIndex(outside, vIdxOutside, dim);
707 DUNE_THROW(Dune::NotImplemented,
"Periodic boundaries for box method for parallel simulations!");
710 const FeCache feCache_;
715 std::size_t numScvf_;
716 std::size_t numBoundaryScvf_;
719 std::vector<bool> boundaryDofIndices_;
722 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
Base class for grid geometries.
Implementation of a boundary face related to primary grid elements (dune intersections).
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices for constant grids.
Definition basegridgeometry.hh:112
void setPeriodic(bool value=true)
Set the periodicity of the grid geometry.
Definition basegridgeometry.hh:169
const GlobalCoordinate & bBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition basegridgeometry.hh:156
Element element(GridIndexType eIdx) const
Get an element from a global element index.
Definition basegridgeometry.hh:142
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices for constant grids.
Definition basegridgeometry.hh:106
const GridView & gridView() const
Return the gridView this grid geometry object lives on.
Definition basegridgeometry.hh:100
void update(const GridView &gridView)
Update all fvElementGeometries (call this after grid adaption).
Definition basegridgeometry.hh:88
const GlobalCoordinate & bBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition basegridgeometry.hh:149
bool isPeriodic() const
Returns if the grid geometry is periodic (at all).
Definition basegridgeometry.hh:162
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition discretization/box/fvelementgeometry.hh:46
std::size_t numBoundaryScvf() const
Definition discretization/box/fvgridgeometry.hh:555
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition discretization/box/fvgridgeometry.hh:598
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition discretization/box/fvgridgeometry.hh:586
GV GridView
export the grid view type
Definition discretization/box/fvgridgeometry.hh:518
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition discretization/box/fvgridgeometry.hh:504
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:508
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition discretization/box/fvgridgeometry.hh:536
static constexpr DiscretizationMethod discMethod
Definition discretization/box/fvgridgeometry.hh:499
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:506
Experimental::BoundaryFace< GV > BoundaryFace
export the boundary face type
Definition discretization/box/fvgridgeometry.hh:510
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvgridgeometry.hh:522
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:590
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition discretization/box/fvgridgeometry.hh:578
BoxGridGeometryCache Cache
Definition discretization/box/fvgridgeometry.hh:624
std::size_t numScvf() const
The total number of sun control volume faces.
Definition discretization/box/fvgridgeometry.hh:550
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition discretization/box/fvgridgeometry.hh:514
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition discretization/box/fvgridgeometry.hh:502
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition discretization/box/fvgridgeometry.hh:582
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvgridgeometry.hh:524
const DofMapper & dofMapper() const
Definition discretization/box/fvgridgeometry.hh:542
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition discretization/box/fvgridgeometry.hh:498
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition discretization/box/fvgridgeometry.hh:520
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:571
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/box/fvgridgeometry.hh:512
std::size_t numDofs() const
The total number of degrees of freedom.
Definition discretization/box/fvgridgeometry.hh:559
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:527
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvgridgeometry.hh:546
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:564
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition discretization/box/fvgridgeometry.hh:516
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition discretization/box/fvgridgeometry.hh:594
GV GridView
export the grid view type
Definition discretization/box/fvgridgeometry.hh:127
Experimental::BoundaryFace< GV > BoundaryFace
export the boundary face type
Definition discretization/box/fvgridgeometry.hh:119
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition discretization/box/fvgridgeometry.hh:207
BoxGridGeometryCache Cache
Definition discretization/box/fvgridgeometry.hh:278
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:173
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:199
std::size_t numDofs() const
The total number of degrees of freedom.
Definition discretization/box/fvgridgeometry.hh:168
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvgridgeometry.hh:155
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition discretization/box/fvgridgeometry.hh:111
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvgridgeometry.hh:131
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition discretization/box/fvgridgeometry.hh:187
std::size_t numBoundaryScvf() const
Definition discretization/box/fvgridgeometry.hh:164
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition discretization/box/fvgridgeometry.hh:107
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvgridgeometry.hh:133
static constexpr DiscretizationMethod discMethod
Definition discretization/box/fvgridgeometry.hh:108
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition discretization/box/fvgridgeometry.hh:125
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:117
const DofMapper & dofMapper() const
Definition discretization/box/fvgridgeometry.hh:151
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition discretization/box/fvgridgeometry.hh:203
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition discretization/box/fvgridgeometry.hh:145
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:136
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/box/fvgridgeometry.hh:121
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition discretization/box/fvgridgeometry.hh:129
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition discretization/box/fvgridgeometry.hh:113
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition discretization/box/fvgridgeometry.hh:191
std::size_t numScvf() const
The total number of sun control volume faces.
Definition discretization/box/fvgridgeometry.hh:159
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition discretization/box/fvgridgeometry.hh:195
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:180
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:115
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition discretization/box/fvgridgeometry.hh:123
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition discretization/box/fvgridgeometry.hh:83
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:60
the sub control volume for the box scheme
Definition discretization/box/subcontrolvolume.hh:58
Class for a boundary face related to primary grid elements (dune intersections).
Definition boundaryface.hh:69
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.
CVFE::DefaultQuadratureTraits< GridView, ScvRule, ScvfRule > BoxQuadratureTraits
Quadrature rule traits for Box discretization.
Definition discretization/box/fvgridgeometry.hh:54
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:38
BaseGridGeometry(std::shared_ptr< BaseImplementation > impl)
Constructor from a BaseImplementation.
Definition basegridgeometry.hh:72
The available discretization methods in Dumux.
Definition cvfelocalresidual.hh:25
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:102
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
Grid properties related to periodicity.
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition discretization/box/fvgridgeometry.hh:65
BoxFVElementGeometry< GridGeometry, enableCache > LocalView
Definition discretization/box/fvgridgeometry.hh:70
BoxSubControlVolume< GridView > SubControlVolume
Definition discretization/box/fvgridgeometry.hh:66
BoxSubControlVolumeFace< GridView > SubControlVolumeFace
Definition discretization/box/fvgridgeometry.hh:67
Quadrature rule traits for discretization schemes.
Definition quadraturerules.hh:85
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:27
unsigned int LocalIndex
Definition indextraits.hh:28
Definition periodicgridtraits.hh:24
Definition periodicgridtraits.hh:23