26#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
31#include <unordered_map>
35#include <dune/geometry/type.hh>
36#include <dune/localfunctions/lagrange/pqkfactory.hh>
52template<
class GG,
bool enableGr
idGeometryCache>
59 using GridView =
typename GG::GridView;
60 static constexpr int dim = GridView::dimension;
61 static constexpr int dimWorld = GridView::dimensionworld;
64 using CoordScalar =
typename GridView::ctype;
65 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
66 using GGCache =
typename GG::Cache;
69 typename GG::SubControlVolume,
70 typename GG::SubControlVolumeFace>;
73 using Element =
typename GridView::template Codim<0>::Entity;
81 static constexpr std::size_t maxNumElementScvs = (1<<dim);
88 : ggCache_(&ggCache) {}
90 [[deprecated(
"This Constructor is deprecated and will be removed after release 3.6. Always use localView(gridGeometry).")]]
97 return ggCache_->scvs(eIdx_)[scvIdx];
103 return ggCache_->scvfs(eIdx_)[scvfIdx];
111 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
114 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
115 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
116 return Dune::IteratorRange<Iter>(s.begin(), s.end());
124 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
127 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
128 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
129 return Dune::IteratorRange<Iter>(s.begin(), s.end());
135 return gridGeometry().feCache().get(element_->type()).localBasis();
141 return ggCache_->scvs(eIdx_).size();
147 return ggCache_->scvfs(eIdx_).size();
157 this->bindElement(element);
158 return std::move(*
this);
165 { this->bindElement(element); }
174 this->bindElement(element);
175 return std::move(*
this);
185 eIdx_ = gridGeometry().elementMapper().index(element);
190 {
return static_cast<bool>(element_); }
194 {
return *element_; }
198 {
return ggCache_->gridGeometry(); }
202 {
return ggCache_->hasBoundaryScvf(eIdx_); }
208 const auto geo = element().geometry();
209 return { Dune::GeometryTypes::cube(dim),
GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
216 const auto geo = element().geometry();
220 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
221 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
222 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
225 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
229 const GGCache* ggCache_;
232 std::optional<Element> element_;
239 using GridView =
typename GG::GridView;
240 static constexpr int dim = GridView::dimension;
241 static constexpr int dimWorld = GridView::dimensionworld;
244 using CoordScalar =
typename GridView::ctype;
245 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
246 using GGCache =
typename GG::Cache;
249 typename GG::SubControlVolume,
250 typename GG::SubControlVolumeFace>;
253 using Element =
typename GridView::template Codim<0>::Entity;
261 static constexpr std::size_t maxNumElementScvs = (1<<dim);
268 : ggCache_(&ggCache) {}
270 [[deprecated(
"This Constructor is deprecated and will be removed after release 3.6. Always use localView(gridGeometry).")]]
277 return scvs_[scvIdx];
283 return scvfs_[scvfIdx];
291 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
294 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
295 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
303 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
306 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
307 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
313 return gridGeometry().feCache().get(element_->type()).localBasis();
325 return scvfs_.size();
335 this->bindElement(element);
336 return std::move(*
this);
343 { this->bindElement(element); }
352 this->bindElement(element);
353 return std::move(*
this);
362 eIdx_ = gridGeometry().elementMapper().index(element);
363 makeElementGeometries_();
368 {
return static_cast<bool>(element_); }
372 {
return *element_; }
376 {
return ggCache_->gridGeometry(); }
380 {
return hasBoundaryScvf_; }
386 const auto geo = element().geometry();
387 return { Dune::GeometryTypes::cube(dim),
GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
394 const auto geo = element().geometry();
398 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
399 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
400 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
403 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
407 void makeElementGeometries_()
409 hasBoundaryScvf_ =
false;
412 const auto& element = *element_;
413 const auto elementGeometry = element.geometry();
414 const auto refElement = referenceElement(elementGeometry);
417 GeometryHelper geometryHelper(elementGeometry);
420 scvs_.resize(elementGeometry.corners());
421 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
424 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
427 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
434 const auto numInnerScvf = geometryHelper.numInteriorScvf();
435 scvfs_.resize(numInnerScvf);
436 scvfBoundaryGeometryKeys_.clear();
438 LocalIndexType scvfLocalIdx = 0;
439 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
442 std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
443 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
445 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
449 std::move(localScvIndices),
454 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
456 if (intersection.boundary() && !intersection.neighbor())
458 const auto isGeometry = intersection.geometry();
459 hasBoundaryScvf_ =
true;
461 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
464 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
465 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
467 scvfs_.emplace_back(geometryHelper,
472 std::move(localScvIndices),
475 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
476 static_cast<LocalIndexType
>(intersection.indexInInside()),
477 static_cast<LocalIndexType
>(isScvfLocalIdx)
489 std::optional<Element> element_;
492 const GGCache* ggCache_;
495 std::vector<SubControlVolume> scvs_;
496 std::vector<SubControlVolumeFace> scvfs_;
497 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
499 bool hasBoundaryScvf_ =
false;
Defines the index types used for grid and local indices.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Class providing iterators over sub control volumes and sub control volume faces of an element.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:53
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:189
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:133
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:172
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:193
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:205
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:112
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:139
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:95
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:145
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:75
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:197
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:201
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:125
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:164
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:213
BoxFVElementGeometry(const GridGeometry &gg)
Definition: discretization/box/fvelementgeometry.hh:91
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:77
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:181
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:101
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:87
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:155
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:73
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:79
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:367
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:323
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:259
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:391
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:379
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:292
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:317
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:359
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:281
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:267
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:383
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:375
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:253
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:304
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:311
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:371
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:275
BoxFVElementGeometry(const GridGeometry &gg)
Definition: discretization/box/fvelementgeometry.hh:271
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:333
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:257
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:350
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:255
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:342