version 3.11-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-FileCopyrightText: 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
30
31namespace Dumux {
32
41template<class GG, bool enableGridGeometryCache>
43
45template<class GG>
46class BoxFVElementGeometry<GG, true>
47{
48 using GridView = typename GG::GridView;
49 static constexpr int dim = GridView::dimension;
50 static constexpr int dimWorld = GridView::dimensionworld;
51 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
52 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
53 using CoordScalar = typename GridView::ctype;
54 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
55 using GGCache = typename GG::Cache;
56 using GeometryHelper = typename GGCache::GeometryHelper;
57
59 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
60 LocalIndexType>;
61public:
63 using Element = typename GridView::template Codim<0>::Entity;
65 using SubControlVolume = typename GG::SubControlVolume;
67 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
69 using GridGeometry = GG;
71 static constexpr std::size_t maxNumElementScvs = (1<<dim);
72
77 BoxFVElementGeometry(const GGCache& ggCache)
78 : ggCache_(&ggCache) {}
79
81 const SubControlVolume& scv(LocalIndexType scvIdx) const
82 {
83 return ggCache_->scvs(eIdx_)[scvIdx];
84 }
85
87 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
88 {
89 return ggCache_->scvfs(eIdx_)[scvfIdx];
90 }
91
97 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
98 scvs(const BoxFVElementGeometry& fvGeometry)
99 {
100 using Iter = typename std::vector<SubControlVolume>::const_iterator;
101 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
102 return Dune::IteratorRange<Iter>(s.begin(), s.end());
103 }
104
110 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
111 scvfs(const BoxFVElementGeometry& fvGeometry)
112 {
113 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
114 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
115 return Dune::IteratorRange<Iter>(s.begin(), s.end());
116 }
117
119 const FeLocalBasis& feLocalBasis() const
120 {
121 return gridGeometry().feCache().get(element_->type()).localBasis();
122 }
123
125 std::size_t numLocalDofs() const
126 {
127 return numScv();
128 }
129
131 std::size_t numScv() const
132 {
133 return ggCache_->scvs(eIdx_).size();
134 }
135
137 std::size_t numScvf() const
138 {
139 return ggCache_->scvfs(eIdx_).size();
140 }
141
148 {
149 this->bindElement(element);
150 return std::move(*this);
151 }
152
156 void bind(const Element& element) &
157 { this->bindElement(element); }
158
165 {
166 this->bindElement(element);
167 return std::move(*this);
168 }
169
173 void bindElement(const Element& element) &
174 {
175 element_ = element;
176 elementGeometry_.emplace(element.geometry());
177 // cache element index
178 eIdx_ = gridGeometry().elementMapper().index(element);
179 }
180
182 bool isBound() const
183 { return static_cast<bool>(element_); }
184
186 const Element& element() const
187 { return *element_; }
188
190 const typename Element::Geometry& elementGeometry() const
191 { return *elementGeometry_; }
192
194 GridIndexType elementIndex() const
195 { return eIdx_; }
196
199 { return ggCache_->gridGeometry(); }
200
202 bool hasBoundaryScvf() const
203 { return ggCache_->hasBoundaryScvf(eIdx_); }
204
206 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
207 {
208 assert(isBound());
209 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
210 }
211
213 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
214 {
215 assert(isBound());
216 const GeometryHelper geometryHelper(*elementGeometry_);
217 if (scvf.boundary())
218 {
219 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
220 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
221 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
222 }
223 else
224 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
225 }
226
228 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
229 {
230 const auto type = fvGeometry.element().type();
231 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
232
233 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
234 }
235
237 template<class LocalDof>
238 friend inline IpData ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
239 {
240 const auto type = fvGeometry.element().type();
241 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
242 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
243
244 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
245 }
246
248 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
249 {
250 // Create ipData that does not automatically calculate the local position but only if it is called
252 [&] (const typename Element::Geometry::GlobalCoordinate& pos)
253 { return fvGeometry.elementGeometry().local(pos); },
254 globalPos
255 };
256 }
257
258private:
259 const GGCache* ggCache_;
260 GridIndexType eIdx_;
261
262 std::optional<Element> element_;
263 std::optional<typename Element::Geometry> elementGeometry_;
264};
265
267template<class GG>
268class BoxFVElementGeometry<GG, false>
269{
270 using GridView = typename GG::GridView;
271 static constexpr int dim = GridView::dimension;
272 static constexpr int dimWorld = GridView::dimensionworld;
273 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
274 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
275 using CoordScalar = typename GridView::ctype;
276 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
277 using GGCache = typename GG::Cache;
278 using GeometryHelper = typename GGCache::GeometryHelper;
280 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
281 LocalIndexType>;
282public:
284 using Element = typename GridView::template Codim<0>::Entity;
286 using SubControlVolume = typename GG::SubControlVolume;
288 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
290 using GridGeometry = GG;
292 static constexpr std::size_t maxNumElementScvs = (1<<dim);
293
298 BoxFVElementGeometry(const GGCache& ggCache)
299 : ggCache_(&ggCache) {}
300
302 const SubControlVolume& scv(LocalIndexType scvIdx) const
303 {
304 return scvs_[scvIdx];
305 }
306
308 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
309 {
310 return scvfs_[scvfIdx];
311 }
312
318 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
319 scvs(const BoxFVElementGeometry& fvGeometry)
320 {
321 using Iter = typename std::vector<SubControlVolume>::const_iterator;
322 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
323 }
324
330 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
331 scvfs(const BoxFVElementGeometry& fvGeometry)
332 {
333 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
334 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
335 }
336
338 const FeLocalBasis& feLocalBasis() const
339 {
340 return gridGeometry().feCache().get(element_->type()).localBasis();
341 }
342
344 std::size_t numLocalDofs() const
345 {
346 return numScv();
347 }
348
350 std::size_t numScv() const
351 {
352 return scvs_.size();
353 }
354
356 std::size_t numScvf() const
357 {
358 return scvfs_.size();
359 }
360
367 {
368 this->bindElement(element);
369 return std::move(*this);
370 }
371
375 void bind(const Element& element) &
376 { this->bindElement(element); }
377
384 {
385 this->bindElement(element);
386 return std::move(*this);
387 }
388
392 void bindElement(const Element& element) &
393 {
394 element_ = element;
395 eIdx_ = gridGeometry().elementMapper().index(element);
396 elementGeometry_.emplace(element.geometry());
397 makeElementGeometries_();
398 }
399
401 bool isBound() const
402 { return static_cast<bool>(element_); }
403
405 const Element& element() const
406 { return *element_; }
407
409 const typename Element::Geometry& elementGeometry() const
410 { return *elementGeometry_; }
411
413 GridIndexType elementIndex() const
414 { return eIdx_; }
415
418 { return ggCache_->gridGeometry(); }
419
421 bool hasBoundaryScvf() const
422 { return hasBoundaryScvf_; }
423
425 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
426 {
427 assert(isBound());
428 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
429 }
430
432 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
433 {
434 assert(isBound());
435 const GeometryHelper geometryHelper(*elementGeometry_);
436 if (scvf.boundary())
437 {
438 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
439 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
440 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
441 }
442 else
443 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
444 }
445
447 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
448 {
449 const auto type = fvGeometry.element().type();
450 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
451
452 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
453 }
454
456 template<class LocalDof>
457 friend inline IpData ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
458 {
459 const auto type = fvGeometry.element().type();
460 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
461 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
462
463 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
464 }
465
467 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
468 {
469 // Create ipData that does not automatically calculate the local position but only if it is called
471 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
472 globalPos
473 };
474 }
475
476private:
477 void makeElementGeometries_()
478 {
479 hasBoundaryScvf_ = false;
480
481 // get the element geometry
482 const auto& element = *element_;
483 const auto& elementGeometry = *elementGeometry_;
484 const auto refElement = referenceElement(elementGeometry);
485
486 // get the sub control volume geometries of this element
487 GeometryHelper geometryHelper(elementGeometry);
488
489 // construct the sub control volumes
490 scvs_.resize(elementGeometry.corners());
491 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
492 {
493 // get associated dof index
494 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
495
496 // add scv to the local container
497 scvs_[scvLocalIdx] = SubControlVolume(
498 geometryHelper.getScvCorners(scvLocalIdx),
499 scvLocalIdx,
500 eIdx_,
501 dofIdxGlobal
502 );
503 }
504
505 // construct the sub control volume faces
506 const auto numInnerScvf = geometryHelper.numInteriorScvf();
507 scvfs_.resize(numInnerScvf);
508 scvfBoundaryGeometryKeys_.clear();
509
510 LocalIndexType scvfLocalIdx = 0;
511 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
512 {
513 // find the local scv indices this scvf is connected to
514 std::array<LocalIndexType, 2> localScvIndices{{
515 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
516 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
517 }};
518
519 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
520 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
521 corners,
522 geometryHelper.normal(corners, localScvIndices),
523 element,
524 scvfLocalIdx,
525 std::move(localScvIndices)
526 );
527 }
528
529 // construct the sub control volume faces on the domain boundary
530 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
531 {
532 if (intersection.boundary() && !intersection.neighbor())
533 {
534 const auto isGeometry = intersection.geometry();
535 hasBoundaryScvf_ = true;
536
537 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
538 {
539 // find the scv this scvf is connected to
540 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
541 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
542
543 scvfs_.emplace_back(
544 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
545 intersection.centerUnitOuterNormal(),
546 intersection,
547 isScvfLocalIdx,
548 scvfLocalIdx,
549 std::move(localScvIndices)
550 );
551
552 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
553 static_cast<LocalIndexType>(intersection.indexInInside()),
554 static_cast<LocalIndexType>(isScvfLocalIdx)
555 }});
556
557 // increment local counter
558 scvfLocalIdx++;
559 }
560 }
561 }
562 }
563
565 GridIndexType eIdx_;
566 std::optional<Element> element_;
567 std::optional<typename Element::Geometry> elementGeometry_;
568
570 const GGCache* ggCache_;
571
573 std::vector<SubControlVolume> scvs_;
574 std::vector<SubControlVolumeFace> scvfs_;
575 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
576
577 bool hasBoundaryScvf_ = false;
578};
579
580} // end namespace Dumux
581
582#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:401
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:356
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:290
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:432
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:421
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:319
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:350
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:413
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:392
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:308
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:298
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:425
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:417
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:284
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:331
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:447
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:338
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:405
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:302
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:467
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:366
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:344
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:457
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:288
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:383
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:409
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:286
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:375
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:182
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:119
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:164
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:186
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:206
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:98
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:131
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:81
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:125
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:137
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:65
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:198
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:202
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:111
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:228
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:190
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:156
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:213
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:248
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:194
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:238
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:67
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:173
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:87
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:77
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:147
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:63
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:69
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:42
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:72
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:50
Classes representing interpolation point data for control-volume finite element schemes.
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