version 3.10-dev
discretization/box/fvelementgeometry.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <utility>
19#include <unordered_map>
20#include <array>
21#include <vector>
22
23#include <dune/geometry/type.hh>
24#include <dune/localfunctions/lagrange/pqkfactory.hh>
25
29
30namespace Dumux {
31
40template<class GG, bool enableGridGeometryCache>
42
44template<class GG>
45class BoxFVElementGeometry<GG, true>
46{
47 using GridView = typename GG::GridView;
48 static constexpr int dim = GridView::dimension;
49 static constexpr int dimWorld = GridView::dimensionworld;
50 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
51 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
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;
56public:
58 using Element = typename GridView::template Codim<0>::Entity;
60 using SubControlVolume = typename GG::SubControlVolume;
62 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
64 using GridGeometry = GG;
66 static constexpr std::size_t maxNumElementScvs = (1<<dim);
67
72 BoxFVElementGeometry(const GGCache& ggCache)
73 : ggCache_(&ggCache) {}
74
76 const SubControlVolume& scv(LocalIndexType scvIdx) const
77 {
78 return ggCache_->scvs(eIdx_)[scvIdx];
79 }
80
82 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
83 {
84 return ggCache_->scvfs(eIdx_)[scvfIdx];
85 }
86
92 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
93 scvs(const BoxFVElementGeometry& fvGeometry)
94 {
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());
98 }
99
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
106 scvfs(const BoxFVElementGeometry& fvGeometry)
107 {
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());
111 }
112
114 const FeLocalBasis& feLocalBasis() const
115 {
116 return gridGeometry().feCache().get(element_->type()).localBasis();
117 }
118
120 std::size_t numScv() const
121 {
122 return ggCache_->scvs(eIdx_).size();
123 }
124
126 std::size_t numScvf() const
127 {
128 return ggCache_->scvfs(eIdx_).size();
129 }
130
137 {
138 this->bindElement(element);
139 return std::move(*this);
140 }
141
145 void bind(const Element& element) &
146 { this->bindElement(element); }
147
154 {
155 this->bindElement(element);
156 return std::move(*this);
157 }
158
162 void bindElement(const Element& element) &
163 {
164 element_ = element;
165 // cache element index
166 eIdx_ = gridGeometry().elementMapper().index(element);
167 }
168
170 bool isBound() const
171 { return static_cast<bool>(element_); }
172
174 const Element& element() const
175 { return *element_; }
176
178 GridIndexType elementIndex() const
179 { return eIdx_; }
180
183 { return ggCache_->gridGeometry(); }
184
186 bool hasBoundaryScvf() const
187 { return ggCache_->hasBoundaryScvf(eIdx_); }
188
190 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
191 {
192 assert(isBound());
193 const auto geo = element().geometry();
194 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
195 }
196
198 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
199 {
200 assert(isBound());
201 const auto geo = element().geometry();
202 const GeometryHelper geometryHelper(geo);
203 if (scvf.boundary())
204 {
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]) };
208 }
209 else
210 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
211 }
212
213private:
214 const GGCache* ggCache_;
215 GridIndexType eIdx_;
216
217 std::optional<Element> element_;
218};
219
221template<class GG>
222class BoxFVElementGeometry<GG, false>
223{
224 using GridView = typename GG::GridView;
225 static constexpr int dim = GridView::dimension;
226 static constexpr int dimWorld = GridView::dimensionworld;
227 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
228 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
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;
233public:
235 using Element = typename GridView::template Codim<0>::Entity;
237 using SubControlVolume = typename GG::SubControlVolume;
239 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
241 using GridGeometry = GG;
243 static constexpr std::size_t maxNumElementScvs = (1<<dim);
244
249 BoxFVElementGeometry(const GGCache& ggCache)
250 : ggCache_(&ggCache) {}
251
253 const SubControlVolume& scv(LocalIndexType scvIdx) const
254 {
255 return scvs_[scvIdx];
256 }
257
259 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
260 {
261 return scvfs_[scvfIdx];
262 }
263
269 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
270 scvs(const BoxFVElementGeometry& fvGeometry)
271 {
272 using Iter = typename std::vector<SubControlVolume>::const_iterator;
273 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
274 }
275
281 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
282 scvfs(const BoxFVElementGeometry& fvGeometry)
283 {
284 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
285 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
286 }
287
289 const FeLocalBasis& feLocalBasis() const
290 {
291 return gridGeometry().feCache().get(element_->type()).localBasis();
292 }
293
295 std::size_t numScv() const
296 {
297 return scvs_.size();
298 }
299
301 std::size_t numScvf() const
302 {
303 return scvfs_.size();
304 }
305
312 {
313 this->bindElement(element);
314 return std::move(*this);
315 }
316
320 void bind(const Element& element) &
321 { this->bindElement(element); }
322
329 {
330 this->bindElement(element);
331 return std::move(*this);
332 }
333
337 void bindElement(const Element& element) &
338 {
339 element_ = element;
340 eIdx_ = gridGeometry().elementMapper().index(element);
341 makeElementGeometries_();
342 }
343
345 bool isBound() const
346 { return static_cast<bool>(element_); }
347
349 const Element& element() const
350 { return *element_; }
351
353 GridIndexType elementIndex() const
354 { return eIdx_; }
355
358 { return ggCache_->gridGeometry(); }
359
361 bool hasBoundaryScvf() const
362 { return hasBoundaryScvf_; }
363
365 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
366 {
367 assert(isBound());
368 const auto geo = element().geometry();
369 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
370 }
371
373 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
374 {
375 assert(isBound());
376 const auto geo = element().geometry();
377 const GeometryHelper geometryHelper(geo);
378 if (scvf.boundary())
379 {
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]) };
383 }
384 else
385 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
386 }
387
388private:
389 void makeElementGeometries_()
390 {
391 hasBoundaryScvf_ = false;
392
393 // get the element geometry
394 const auto& element = *element_;
395 const auto elementGeometry = element.geometry();
396 const auto refElement = referenceElement(elementGeometry);
397
398 // get the sub control volume geometries of this element
399 GeometryHelper geometryHelper(elementGeometry);
400
401 // construct the sub control volumes
402 scvs_.resize(elementGeometry.corners());
403 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
404 {
405 // get associated dof index
406 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
407
408 // add scv to the local container
409 scvs_[scvLocalIdx] = SubControlVolume(
410 geometryHelper.getScvCorners(scvLocalIdx),
411 scvLocalIdx,
412 eIdx_,
413 dofIdxGlobal
414 );
415 }
416
417 // construct the sub control volume faces
418 const auto numInnerScvf = geometryHelper.numInteriorScvf();
419 scvfs_.resize(numInnerScvf);
420 scvfBoundaryGeometryKeys_.clear();
421
422 LocalIndexType scvfLocalIdx = 0;
423 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
424 {
425 // find the local scv indices this scvf is connected to
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))});
428
429 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
430 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
431 corners,
432 geometryHelper.normal(corners, localScvIndices),
433 element,
434 elementGeometry,
435 scvfLocalIdx,
436 std::move(localScvIndices),
437 false
438 );
439 }
440
441 // construct the sub control volume faces on the domain boundary
442 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
443 {
444 if (intersection.boundary() && !intersection.neighbor())
445 {
446 const auto isGeometry = intersection.geometry();
447 hasBoundaryScvf_ = true;
448
449 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
450 {
451 // find the scv this scvf is connected to
452 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
453 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
454
455 scvfs_.emplace_back(
456 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
457 intersection.centerUnitOuterNormal(),
458 intersection,
459 isGeometry,
460 isScvfLocalIdx,
461 scvfLocalIdx,
462 std::move(localScvIndices),
463 true
464 );
465
466 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
467 static_cast<LocalIndexType>(intersection.indexInInside()),
468 static_cast<LocalIndexType>(isScvfLocalIdx)
469 }});
470
471 // increment local counter
472 scvfLocalIdx++;
473 }
474 }
475 }
476 }
477
479 GridIndexType eIdx_;
480 std::optional<Element> element_;
481
483 const GGCache* ggCache_;
484
486 std::vector<SubControlVolume> scvs_;
487 std::vector<SubControlVolumeFace> scvfs_;
488 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
489
490 bool hasBoundaryScvf_ = false;
491};
492
493} // end namespace Dumux
494
495#endif
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
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
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