14#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
19#include <unordered_map>
23#include <dune/geometry/type.hh>
24#include <dune/localfunctions/lagrange/pqkfactory.hh>
40template<
class GG,
bool enableGr
idGeometryCache>
47 using GridView =
typename GG::GridView;
48 static constexpr int dim = GridView::dimension;
49 static constexpr int dimWorld = GridView::dimensionworld;
52 using CoordScalar =
typename GridView::ctype;
53 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
54 using GGCache =
typename GG::Cache;
55 using GeometryHelper =
typename GGCache::GeometryHelper;
58 using Element =
typename GridView::template Codim<0>::Entity;
66 static constexpr std::size_t maxNumElementScvs = (1<<dim);
73 : ggCache_(&ggCache) {}
78 return ggCache_->scvs(eIdx_)[scvIdx];
84 return ggCache_->scvfs(eIdx_)[scvfIdx];
92 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
95 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
96 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
97 return Dune::IteratorRange<Iter>(s.begin(), s.end());
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
108 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
110 return Dune::IteratorRange<Iter>(s.begin(), s.end());
116 return gridGeometry().feCache().get(element_->type()).localBasis();
122 return ggCache_->scvs(eIdx_).size();
128 return ggCache_->scvfs(eIdx_).size();
138 this->bindElement(element);
139 return std::move(*
this);
146 { this->bindElement(element); }
155 this->bindElement(element);
156 return std::move(*
this);
166 eIdx_ = gridGeometry().elementMapper().index(element);
171 {
return static_cast<bool>(element_); }
175 {
return *element_; }
183 {
return ggCache_->gridGeometry(); }
187 {
return ggCache_->hasBoundaryScvf(eIdx_); }
193 const auto geo = element().geometry();
194 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
201 const auto geo = element().geometry();
202 const GeometryHelper geometryHelper(geo);
205 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
206 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
207 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
210 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
214 const GGCache* ggCache_;
217 std::optional<Element> element_;
224 using GridView =
typename GG::GridView;
225 static constexpr int dim = GridView::dimension;
226 static constexpr int dimWorld = GridView::dimensionworld;
229 using CoordScalar =
typename GridView::ctype;
230 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
231 using GGCache =
typename GG::Cache;
232 using GeometryHelper =
typename GGCache::GeometryHelper;
235 using Element =
typename GridView::template Codim<0>::Entity;
243 static constexpr std::size_t maxNumElementScvs = (1<<dim);
250 : ggCache_(&ggCache) {}
255 return scvs_[scvIdx];
261 return scvfs_[scvfIdx];
269 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
272 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
273 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
281 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
284 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
285 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
291 return gridGeometry().feCache().get(element_->type()).localBasis();
303 return scvfs_.size();
313 this->bindElement(element);
314 return std::move(*
this);
321 { this->bindElement(element); }
330 this->bindElement(element);
331 return std::move(*
this);
340 eIdx_ = gridGeometry().elementMapper().index(element);
341 makeElementGeometries_();
346 {
return static_cast<bool>(element_); }
350 {
return *element_; }
358 {
return ggCache_->gridGeometry(); }
362 {
return hasBoundaryScvf_; }
368 const auto geo = element().geometry();
369 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
376 const auto geo = element().geometry();
377 const GeometryHelper geometryHelper(geo);
380 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
381 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
382 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
385 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
389 void makeElementGeometries_()
391 hasBoundaryScvf_ =
false;
394 const auto& element = *element_;
395 const auto elementGeometry = element.geometry();
396 const auto refElement = referenceElement(elementGeometry);
399 GeometryHelper geometryHelper(elementGeometry);
402 scvs_.resize(elementGeometry.corners());
403 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
406 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
409 scvs_[scvLocalIdx] = SubControlVolume(
410 geometryHelper.getScvCorners(scvLocalIdx),
418 const auto numInnerScvf = geometryHelper.numInteriorScvf();
419 scvfs_.resize(numInnerScvf);
420 scvfBoundaryGeometryKeys_.clear();
422 LocalIndexType scvfLocalIdx = 0;
423 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
426 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
427 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
429 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
430 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
432 geometryHelper.normal(corners, localScvIndices),
436 std::move(localScvIndices),
442 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
444 if (intersection.boundary() && !intersection.neighbor())
446 const auto isGeometry = intersection.geometry();
447 hasBoundaryScvf_ =
true;
449 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
452 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
453 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
456 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
457 intersection.centerUnitOuterNormal(),
462 std::move(localScvIndices),
466 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
467 static_cast<LocalIndexType
>(intersection.indexInInside()),
468 static_cast<LocalIndexType
>(isScvfLocalIdx)
480 std::optional<Element> element_;
483 const GGCache* ggCache_;
486 std::vector<SubControlVolume> scvs_;
487 std::vector<SubControlVolumeFace> scvfs_;
488 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
490 bool hasBoundaryScvf_ =
false;
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:345
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:301
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:241
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:373
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:361
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:270
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:295
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:353
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:337
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:259
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:249
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:365
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:357
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:235
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:282
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:289
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:349
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:253
BoxFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:311
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:239
BoxFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:328
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:237
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:320
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:170
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:114
BoxFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:153
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:174
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:190
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:93
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:120
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:76
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:126
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:60
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:182
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:186
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:106
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:145
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:198
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:178
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:62
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:162
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:82
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:72
BoxFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/box/fvelementgeometry.hh:136
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:58
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:64
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:41
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28