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
31
32namespace Dumux {
33
42template<class GG, bool enableGridGeometryCache>
44
46template<class GG>
47class BoxFVElementGeometry<GG, true>
48{
49 using GridView = typename GG::GridView;
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
52 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
53 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
54 using CoordScalar = typename GridView::ctype;
55 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
56 using GGCache = typename GG::Cache;
57 using GeometryHelper = typename GGCache::GeometryHelper;
58
60 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
61 LocalIndexType>;
62public:
64 using Element = typename GridView::template Codim<0>::Entity;
66 using SubControlVolume = typename GG::SubControlVolume;
68 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
70 using GridGeometry = GG;
72 static constexpr std::size_t maxNumElementScvs = (1<<dim);
73
78 BoxFVElementGeometry(const GGCache& ggCache)
79 : ggCache_(&ggCache) {}
80
82 const SubControlVolume& scv(LocalIndexType scvIdx) const
83 {
84 return ggCache_->scvs(eIdx_)[scvIdx];
85 }
86
88 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
89 {
90 return ggCache_->scvfs(eIdx_)[scvfIdx];
91 }
92
98 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
99 scvs(const BoxFVElementGeometry& fvGeometry)
100 {
101 using Iter = typename std::vector<SubControlVolume>::const_iterator;
102 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
103 return Dune::IteratorRange<Iter>(s.begin(), s.end());
104 }
105
111 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
112 scvfs(const BoxFVElementGeometry& fvGeometry)
113 {
114 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
115 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
116 return Dune::IteratorRange<Iter>(s.begin(), s.end());
117 }
118
120 template<class Intersection>
121 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const Intersection& intersection)
122 {
123 return Dune::transformedRangeView(
124 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
125 [&](const auto i)
126 {
127 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
128 return CVFE::LocalDof
129 {
130 static_cast<LocalIndexType>(localDofIdx),
131 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
132 static_cast<GridIndexType>(fvGeometry.elementIndex())
133 };
134 }
135 );
136 }
137
139 const FeLocalBasis& feLocalBasis() const
140 {
141 return gridGeometry().feCache().get(element_->type()).localBasis();
142 }
143
145 std::size_t numLocalDofs() const
146 {
147 return numScv();
148 }
149
151 std::size_t numScv() const
152 {
153 return ggCache_->scvs(eIdx_).size();
154 }
155
157 std::size_t numScvf() const
158 {
159 return ggCache_->scvfs(eIdx_).size();
160 }
161
168 {
169 this->bindElement(element);
170 return std::move(*this);
171 }
172
176 void bind(const Element& element) &
177 { this->bindElement(element); }
178
185 {
186 this->bindElement(element);
187 return std::move(*this);
188 }
189
193 void bindElement(const Element& element) &
194 {
195 element_ = element;
196 elementGeometry_.emplace(element.geometry());
197 // cache element index
198 eIdx_ = gridGeometry().elementMapper().index(element);
199 }
200
202 bool isBound() const
203 { return static_cast<bool>(element_); }
204
206 const Element& element() const
207 { return *element_; }
208
210 const typename Element::Geometry& elementGeometry() const
211 { return *elementGeometry_; }
212
214 GridIndexType elementIndex() const
215 { return eIdx_; }
216
218 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
219 {
220 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
221 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
222 }
223
226 { return ggCache_->gridGeometry(); }
227
229 bool hasBoundaryScvf() const
230 { return ggCache_->hasBoundaryScvf(eIdx_); }
231
233 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
234 {
235 assert(isBound());
236 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
237 }
238
240 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
241 {
242 assert(isBound());
243 const GeometryHelper geometryHelper(*elementGeometry_);
244 if (scvf.boundary())
245 {
246 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
247 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
248 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
249 }
250 else
251 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
252 }
253
255 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
256 {
257 const auto type = fvGeometry.element().type();
258 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
259
260 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
261 }
262
264 template<class LocalDof>
265 friend inline IpData ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
266 {
267 const auto type = fvGeometry.element().type();
268 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
269 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
270
271 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
272 }
273
275 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
276 {
277 // Create ipData that does not automatically calculate the local position but only if it is called
279 [&] (const typename Element::Geometry::GlobalCoordinate& pos)
280 { return fvGeometry.elementGeometry().local(pos); },
281 globalPos
282 };
283 }
284
285private:
286 const GGCache* ggCache_;
287 GridIndexType eIdx_;
288
289 std::optional<Element> element_;
290 std::optional<typename Element::Geometry> elementGeometry_;
291};
292
294template<class GG>
295class BoxFVElementGeometry<GG, false>
296{
297 using GridView = typename GG::GridView;
298 static constexpr int dim = GridView::dimension;
299 static constexpr int dimWorld = GridView::dimensionworld;
300 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
301 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
302 using CoordScalar = typename GridView::ctype;
303 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
304 using GGCache = typename GG::Cache;
305 using GeometryHelper = typename GGCache::GeometryHelper;
307 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
308 LocalIndexType>;
309public:
311 using Element = typename GridView::template Codim<0>::Entity;
313 using SubControlVolume = typename GG::SubControlVolume;
315 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
317 using GridGeometry = GG;
319 static constexpr std::size_t maxNumElementScvs = (1<<dim);
320
325 BoxFVElementGeometry(const GGCache& ggCache)
326 : ggCache_(&ggCache) {}
327
329 const SubControlVolume& scv(LocalIndexType scvIdx) const
330 {
331 return scvs_[scvIdx];
332 }
333
335 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
336 {
337 return scvfs_[scvfIdx];
338 }
339
345 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
346 scvs(const BoxFVElementGeometry& fvGeometry)
347 {
348 using Iter = typename std::vector<SubControlVolume>::const_iterator;
349 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
350 }
351
357 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
358 scvfs(const BoxFVElementGeometry& fvGeometry)
359 {
360 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
361 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
362 }
363
365 template<class Intersection>
366 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const Intersection& intersection)
367 {
368 return Dune::transformedRangeView(
369 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
370 [&](const auto i)
371 {
372 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
373 return CVFE::LocalDof
374 {
375 static_cast<LocalIndexType>(localDofIdx),
376 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
377 static_cast<GridIndexType>(fvGeometry.elementIndex())
378 };
379 }
380 );
381 }
382
383
385 const FeLocalBasis& feLocalBasis() const
386 {
387 return gridGeometry().feCache().get(element_->type()).localBasis();
388 }
389
391 std::size_t numLocalDofs() const
392 {
393 return numScv();
394 }
395
397 std::size_t numScv() const
398 {
399 return scvs_.size();
400 }
401
403 std::size_t numScvf() const
404 {
405 return scvfs_.size();
406 }
407
414 {
415 this->bindElement(element);
416 return std::move(*this);
417 }
418
422 void bind(const Element& element) &
423 { this->bindElement(element); }
424
431 {
432 this->bindElement(element);
433 return std::move(*this);
434 }
435
439 void bindElement(const Element& element) &
440 {
441 element_ = element;
442 eIdx_ = gridGeometry().elementMapper().index(element);
443 elementGeometry_.emplace(element.geometry());
444 makeElementGeometries_();
445 }
446
448 bool isBound() const
449 { return static_cast<bool>(element_); }
450
452 const Element& element() const
453 { return *element_; }
454
456 const typename Element::Geometry& elementGeometry() const
457 { return *elementGeometry_; }
458
460 GridIndexType elementIndex() const
461 { return eIdx_; }
462
464 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
465 {
466 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
467 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
468 }
469
472 { return ggCache_->gridGeometry(); }
473
475 bool hasBoundaryScvf() const
476 { return hasBoundaryScvf_; }
477
479 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
480 {
481 assert(isBound());
482 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
483 }
484
486 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
487 {
488 assert(isBound());
489 const GeometryHelper geometryHelper(*elementGeometry_);
490 if (scvf.boundary())
491 {
492 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
493 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
494 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
495 }
496 else
497 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
498 }
499
501 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
502 {
503 const auto type = fvGeometry.element().type();
504 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
505
506 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
507 }
508
510 template<class LocalDof>
511 friend inline IpData ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
512 {
513 const auto type = fvGeometry.element().type();
514 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
515 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
516
517 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
518 }
519
521 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
522 {
523 // Create ipData that does not automatically calculate the local position but only if it is called
525 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
526 globalPos
527 };
528 }
529
530private:
531 void makeElementGeometries_()
532 {
533 hasBoundaryScvf_ = false;
534
535 // get the element geometry
536 const auto& element = *element_;
537 const auto& elementGeometry = *elementGeometry_;
538 const auto refElement = referenceElement(elementGeometry);
539
540 // get the sub control volume geometries of this element
541 GeometryHelper geometryHelper(elementGeometry);
542
543 // construct the sub control volumes
544 scvs_.resize(elementGeometry.corners());
545 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
546 {
547 // get associated dof index
548 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
549
550 // add scv to the local container
551 scvs_[scvLocalIdx] = SubControlVolume(
552 geometryHelper.getScvCorners(scvLocalIdx),
553 scvLocalIdx,
554 eIdx_,
555 dofIdxGlobal
556 );
557 }
558
559 // construct the sub control volume faces
560 const auto numInnerScvf = geometryHelper.numInteriorScvf();
561 scvfs_.resize(numInnerScvf);
562 scvfBoundaryGeometryKeys_.clear();
563
564 LocalIndexType scvfLocalIdx = 0;
565 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
566 {
567 // find the local scv indices this scvf is connected to
568 std::array<LocalIndexType, 2> localScvIndices{{
569 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
570 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
571 }};
572
573 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
574 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
575 corners,
576 geometryHelper.normal(corners, localScvIndices),
577 element,
578 scvfLocalIdx,
579 std::move(localScvIndices)
580 );
581 }
582
583 // construct the sub control volume faces on the domain boundary
584 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
585 {
586 if (intersection.boundary() && !intersection.neighbor())
587 {
588 const auto isGeometry = intersection.geometry();
589 hasBoundaryScvf_ = true;
590
591 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
592 {
593 // find the scv this scvf is connected to
594 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
595 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
596
597 scvfs_.emplace_back(
598 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
599 intersection.centerUnitOuterNormal(),
600 intersection,
601 isScvfLocalIdx,
602 scvfLocalIdx,
603 std::move(localScvIndices)
604 );
605
606 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
607 static_cast<LocalIndexType>(intersection.indexInInside()),
608 static_cast<LocalIndexType>(isScvfLocalIdx)
609 }});
610
611 // increment local counter
612 scvfLocalIdx++;
613 }
614 }
615 }
616 }
617
619 GridIndexType eIdx_;
620 std::optional<Element> element_;
621 std::optional<typename Element::Geometry> elementGeometry_;
622
624 const GGCache* ggCache_;
625
627 std::vector<SubControlVolume> scvs_;
628 std::vector<SubControlVolumeFace> scvfs_;
629 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
630
631 bool hasBoundaryScvf_ = false;
632};
633
634} // end namespace Dumux
635
636#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:448
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:366
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:403
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:317
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:486
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:475
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:346
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:397
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:460
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:439
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:464
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:335
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:325
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:479
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:471
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:311
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:358
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:501
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:385
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:452
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:329
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:521
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:413
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:391
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:511
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:315
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:430
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:456
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:313
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:422
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:202
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:139
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:121
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:184
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:206
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:233
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:99
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:151
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:82
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:145
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:157
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:66
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:225
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:229
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:112
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:255
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:210
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:176
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:218
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:240
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:275
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:214
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:265
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:68
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:193
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:88
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:78
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:167
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:64
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:70
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:43
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.
Class representing dofs on elements for control-volume finite element schemes.
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