version 3.11-dev
Loading...
Searching...
No Matches
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 <span>
19#include <utility>
20#include <unordered_map>
21#include <array>
22#include <vector>
23
24#include <dune/geometry/type.hh>
25#include <dune/localfunctions/lagrange/pqkfactory.hh>
26#include <dune/common/reservedvector.hh>
27
34
35namespace Dumux {
36
45template<class GG, bool enableGridGeometryCache>
47
49template<class GG>
50class BoxFVElementGeometry<GG, true>
51{
52 using GridView = typename GG::GridView;
53 static constexpr int dim = GridView::dimension;
54 static constexpr int dimWorld = GridView::dimensionworld;
55 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
56 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
57 using CoordScalar = typename GridView::ctype;
58 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
59 using GGCache = typename GG::Cache;
60 using GeometryHelper = typename GGCache::GeometryHelper;
61
62 using BaseIpData = CVFE::InterpolationPointData<
63 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
64 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
65 >;
66
67public:
69 using Element = typename GridView::template Codim<0>::Entity;
71 using SubControlVolume = typename GG::SubControlVolume;
73 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
75 using GridGeometry = GG;
77 using BoundaryFace = typename GG::BoundaryFace;
79 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
81 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
83 static constexpr std::size_t maxNumElementScvs = (1<<dim);
84
89 BoxFVElementGeometry(const GGCache& ggCache)
90 : ggCache_(&ggCache) {}
91
93 const SubControlVolume& scv(LocalIndexType scvIdx) const
94 {
95 return ggCache_->scvs(eIdx_)[scvIdx];
96 }
97
99 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
100 {
101 return ggCache_->scvfs(eIdx_)[scvfIdx];
102 }
103
109 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
110 scvs(const BoxFVElementGeometry& fvGeometry)
111 {
112 using Iter = typename std::vector<SubControlVolume>::const_iterator;
113 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
114 return Dune::IteratorRange<Iter>(s.begin(), s.end());
115 }
116
122 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
123 scvfs(const BoxFVElementGeometry& fvGeometry)
124 {
125 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
126 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
127 return Dune::IteratorRange<Iter>(s.begin(), s.end());
128 }
129
132 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
134 {
135 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
136 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
137 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[boundaryFace.index()];
138 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
139 }
140
143 friend inline std::ranges::view auto
145 {
146 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
147 return std::ranges::views::all(v);
148 }
149
151 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const BoundaryFace& boundaryFace)
152 {
153 return Dune::transformedRangeView(
154 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(boundaryFace.intersectionIndex(), 1, dim)),
155 [&](const auto i)
156 {
157 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(boundaryFace.intersectionIndex(), 1, i, dim);
158 return CVFE::LocalDof
159 {
160 static_cast<LocalIndexType>(localDofIdx),
161 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
162 static_cast<GridIndexType>(fvGeometry.elementIndex())
163 };
164 }
165 );
166 }
167
169 const FeLocalBasis& feLocalBasis() const
170 {
171 return gridGeometry().feCache().get(element_->type()).localBasis();
172 }
173
175 std::size_t numLocalDofs() const
176 {
177 return numScv();
178 }
179
181 std::size_t numScv() const
182 {
183 return ggCache_->scvs(eIdx_).size();
184 }
185
187 std::size_t numScvf() const
188 {
189 return ggCache_->scvfs(eIdx_).size();
190 }
191
198 {
199 this->bindElement(element);
200 return std::move(*this);
201 }
202
206 void bind(const Element& element) &
207 { this->bindElement(element); }
208
215 {
216 this->bindElement(element);
217 return std::move(*this);
218 }
219
224 {
225 element_ = element;
226 elementGeometry_.emplace(element.geometry());
227 // cache element index
228 eIdx_ = gridGeometry().elementMapper().index(element);
229 }
230
232 bool isBound() const
233 { return static_cast<bool>(element_); }
234
236 const Element& element() const
237 { return *element_; }
238
240 const typename Element::Geometry& elementGeometry() const
241 { return *elementGeometry_; }
242
244 GridIndexType elementIndex() const
245 { return eIdx_; }
246
249 {
250 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
251 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
252 }
253
256 { return ggCache_->gridGeometry(); }
257
259 bool hasBoundaryScvf() const
260 { return ggCache_->hasBoundaryScvf(eIdx_); }
261
263 bool hasBoundaryFaces() const
264 { return hasBoundaryScvf(); }
265
267 const BoundaryFace& boundaryFace(LocalIndexType bfIdx) const
268 { return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
269
271 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
272 {
273 assert(isBound());
274 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
275 }
276
278 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
279 {
280 assert(isBound());
281 const GeometryHelper geometryHelper(*elementGeometry_);
282 if (scvf.boundary())
283 {
284 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
285 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
286 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
287 }
288 else
289 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
290 }
291
293 typename BoundaryFace::Traits::Geometry geometry(const BoundaryFace& boundaryFace) const
294 {
295 assert(isBound());
296 const auto& elemGeo = elementGeometry();
297 const auto faceGeoInRef = referenceElement(elemGeo).template geometry<1>(boundaryFace.intersectionIndex());
298 typename BoundaryFace::Traits::CornerStorage corners;
299 for (int i = 0; i < faceGeoInRef.corners(); ++i)
300 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
301 return { faceGeoInRef.type(), corners };
302 }
303
305 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
306 {
307 const auto type = fvGeometry.element().type();
308 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
309
310 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
311 }
312
314 template<class LocalDof>
315 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
316 {
317 const auto type = fvGeometry.element().type();
318 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
319 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
320
321 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
322 }
323
325 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
326 {
327 // Create ipData that does not automatically calculate the local position but only if it is called
329 [&] (const typename Element::Geometry::GlobalCoordinate& pos)
330 { return fvGeometry.elementGeometry().local(pos); },
331 globalPos
332 };
333 }
334
336 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
337 {
339 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
340 }
341
342private:
343 const GGCache* ggCache_;
344 GridIndexType eIdx_;
345
346 std::optional<Element> element_;
347 std::optional<typename Element::Geometry> elementGeometry_;
348};
349
351template<class GG>
352class BoxFVElementGeometry<GG, false>
353{
354 using GridView = typename GG::GridView;
355 static constexpr int dim = GridView::dimension;
356 static constexpr int dimWorld = GridView::dimensionworld;
357 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
358 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
359 using CoordScalar = typename GridView::ctype;
360 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
361 using GGCache = typename GG::Cache;
362 using GeometryHelper = typename GGCache::GeometryHelper;
363
364 using BaseIpData = CVFE::InterpolationPointData<
365 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
366 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
367 >;
368
369public:
371 using Element = typename GridView::template Codim<0>::Entity;
373 using SubControlVolume = typename GG::SubControlVolume;
375 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
377 using GridGeometry = GG;
379 using BoundaryFace = typename GG::BoundaryFace;
381 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
383 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
385 static constexpr std::size_t maxNumElementScvs = (1<<dim);
386
391 BoxFVElementGeometry(const GGCache& ggCache)
392 : ggCache_(&ggCache) {}
393
395 const SubControlVolume& scv(LocalIndexType scvIdx) const
396 {
397 return scvs_[scvIdx];
398 }
399
401 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
402 {
403 return scvfs_[scvfIdx];
404 }
405
411 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
412 scvs(const BoxFVElementGeometry& fvGeometry)
413 {
414 using Iter = typename std::vector<SubControlVolume>::const_iterator;
415 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
416 }
417
423 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
424 scvfs(const BoxFVElementGeometry& fvGeometry)
425 {
426 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
427 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
428 }
429
432 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
434 {
435 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
436 const auto& range = fvGeometry.boundaryFaceScvfRanges_[boundaryFace.index()];
437 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin() + range[0], fvGeometry.scvfs_.begin() + range[0] + range[1]);
438 }
439
442 friend inline std::ranges::view auto
444 { return std::ranges::views::all(fvGeometry.boundaryFaces_); }
445
447 friend inline auto localDofs(const BoxFVElementGeometry& fvGeometry, const BoundaryFace& boundaryFace)
448 {
449 return Dune::transformedRangeView(
450 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(boundaryFace.intersectionIndex(), 1, dim)),
451 [&](const auto i)
452 {
453 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(boundaryFace.intersectionIndex(), 1, i, dim);
454 return CVFE::LocalDof
455 {
456 static_cast<LocalIndexType>(localDofIdx),
457 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
458 static_cast<GridIndexType>(fvGeometry.elementIndex())
459 };
460 }
461 );
462 }
463
464
466 const FeLocalBasis& feLocalBasis() const
467 {
468 return gridGeometry().feCache().get(element_->type()).localBasis();
469 }
470
472 std::size_t numLocalDofs() const
473 {
474 return numScv();
475 }
476
478 std::size_t numScv() const
479 {
480 return scvs_.size();
481 }
482
484 std::size_t numScvf() const
485 {
486 return scvfs_.size();
487 }
488
495 {
496 this->bindElement(element);
497 return std::move(*this);
498 }
499
503 void bind(const Element& element) &
504 { this->bindElement(element); }
505
512 {
513 this->bindElement(element);
514 return std::move(*this);
515 }
516
521 {
522 element_ = element;
523 eIdx_ = gridGeometry().elementMapper().index(element);
524 elementGeometry_.emplace(element.geometry());
525 makeElementGeometries_();
526 }
527
529 bool isBound() const
530 { return static_cast<bool>(element_); }
531
533 const Element& element() const
534 { return *element_; }
535
537 const typename Element::Geometry& elementGeometry() const
538 { return *elementGeometry_; }
539
541 GridIndexType elementIndex() const
542 { return eIdx_; }
543
546 {
547 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
548 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
549 }
550
553 { return ggCache_->gridGeometry(); }
554
556 bool hasBoundaryScvf() const
557 { return hasBoundaryScvf_; }
558
560 bool hasBoundaryFaces() const
561 { return hasBoundaryScvf(); }
562
564 const BoundaryFace& boundaryFace(LocalIndexType bfIdx) const
565 { return boundaryFaces_[bfIdx]; }
566
568 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
569 {
570 assert(isBound());
571 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
572 }
573
575 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
576 {
577 assert(isBound());
578 const GeometryHelper geometryHelper(*elementGeometry_);
579 if (scvf.boundary())
580 {
581 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
582 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
583 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
584 }
585 else
586 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
587 }
588
590 typename BoundaryFace::Traits::Geometry geometry(const BoundaryFace& boundaryFace) const
591 {
592 assert(isBound());
593 const auto& elemGeo = elementGeometry();
594 const auto faceGeoInRef = referenceElement(elemGeo).template geometry<1>(boundaryFace.intersectionIndex());
595 typename BoundaryFace::Traits::CornerStorage corners;
596 for (int i = 0; i < faceGeoInRef.corners(); ++i)
597 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
598 return { faceGeoInRef.type(), corners };
599 }
600
602 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolume& scv)
603 {
604 const auto type = fvGeometry.element().type();
605 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
606
607 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
608 }
609
611 template<class LocalDof>
612 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const LocalDof& localDof)
613 {
614 const auto type = fvGeometry.element().type();
615 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
616 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
617
618 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
619 }
620
622 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
623 {
624 // Create ipData that does not automatically calculate the local position but only if it is called
626 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
627 globalPos
628 };
629 }
630
632 friend inline auto ipData(const BoxFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
633 {
635 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
636 }
637
638private:
639 void makeElementGeometries_()
640 {
641 hasBoundaryScvf_ = false;
642 boundaryFaces_.clear();
643 boundaryFaceScvfRanges_.clear();
644
645 // get the element geometry
646 const auto& element = *element_;
647 const auto& elementGeometry = *elementGeometry_;
648 const auto refElement = referenceElement(elementGeometry);
649
650 // get the sub control volume geometries of this element
651 GeometryHelper geometryHelper(elementGeometry);
652
653 // construct the sub control volumes
654 scvs_.resize(elementGeometry.corners());
655 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
656 {
657 // get associated dof index
658 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
659
660 // add scv to the local container
661 scvs_[scvLocalIdx] = SubControlVolume(
662 geometryHelper.getScvCorners(scvLocalIdx),
663 scvLocalIdx,
664 eIdx_,
665 dofIdxGlobal
666 );
667 }
668
669 // construct the sub control volume faces
670 const auto numInnerScvf = geometryHelper.numInteriorScvf();
671 scvfs_.resize(numInnerScvf);
672 scvfBoundaryGeometryKeys_.clear();
673
674 LocalIndexType scvfLocalIdx = 0;
675 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
676 {
677 // find the local scv indices this scvf is connected to
678 std::array<LocalIndexType, 2> localScvIndices{{
679 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
680 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
681 }};
682
683 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
684 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
685 corners,
686 geometryHelper.normal(corners, localScvIndices),
687 element,
688 scvfLocalIdx,
689 std::move(localScvIndices)
690 );
691 }
692
693 // construct the sub control volume faces on the domain boundary
694 LocalIndexType numBoundaryFaces = 0;
695 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
696 {
697 if (intersection.boundary() && !intersection.neighbor())
698 {
699 const auto isGeometry = intersection.geometry();
700 hasBoundaryScvf_ = true;
701
702 // add one boundary face per boundary intersection
703 boundaryFaces_.push_back(BoundaryFace{
704 isGeometry.center(),
705 isGeometry.volume(),
706 intersection.centerUnitOuterNormal(),
707 numBoundaryFaces++,
708 static_cast<LocalIndexType>(intersection.indexInInside()),
709 typename BoundaryFace::Traits::BoundaryFlag{intersection}
710 });
711
712 // record the scvf subrange for this boundary face: {offset, count}
713 boundaryFaceScvfRanges_.push_back(std::array<LocalIndexType, 2>{{
714 scvfLocalIdx,
715 static_cast<LocalIndexType>(isGeometry.corners())
716 }});
717
718 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
719 {
720 // find the scv this scvf is connected to
721 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
722 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
723
724 scvfs_.emplace_back(
725 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
726 intersection.centerUnitOuterNormal(),
727 intersection,
728 isScvfLocalIdx,
729 scvfLocalIdx,
730 std::move(localScvIndices)
731 );
732
733 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
734 static_cast<LocalIndexType>(intersection.indexInInside()),
735 static_cast<LocalIndexType>(isScvfLocalIdx)
736 }});
737
738 // increment local counter
739 scvfLocalIdx++;
740 }
741 }
742 }
743 }
744
746 GridIndexType eIdx_;
747 std::optional<Element> element_;
748 std::optional<typename Element::Geometry> elementGeometry_;
749
751 const GGCache* ggCache_;
752
754 std::vector<SubControlVolume> scvs_;
755 std::vector<SubControlVolumeFace> scvfs_;
756 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
757 Dune::ReservedVector<BoundaryFace, 2*dim> boundaryFaces_;
758 Dune::ReservedVector<std::array<LocalIndexType, 2>, 2*dim> boundaryFaceScvfRanges_;
759
760 bool hasBoundaryScvf_ = false;
761};
762
763} // end namespace Dumux
764
765#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:529
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/box/fvelementgeometry.hh:484
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/box/fvelementgeometry.hh:377
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/box/fvelementgeometry.hh:575
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/box/fvelementgeometry.hh:556
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:412
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvelementgeometry.hh:478
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition discretization/box/fvelementgeometry.hh:541
void bindElement(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:520
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/box/fvelementgeometry.hh:545
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/box/fvelementgeometry.hh:401
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/box/fvelementgeometry.hh:391
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvelementgeometry.hh:383
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/box/fvelementgeometry.hh:568
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/box/fvelementgeometry.hh:552
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/box/fvelementgeometry.hh:433
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/box/fvelementgeometry.hh:371
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:424
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/box/fvelementgeometry.hh:602
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/box/fvelementgeometry.hh:590
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/box/fvelementgeometry.hh:466
const Element & element() const
The bound element.
Definition discretization/box/fvelementgeometry.hh:533
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/box/fvelementgeometry.hh:395
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/box/fvelementgeometry.hh:622
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:494
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/box/fvelementgeometry.hh:632
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/box/fvelementgeometry.hh:447
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/box/fvelementgeometry.hh:612
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/box/fvelementgeometry.hh:381
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/box/fvelementgeometry.hh:472
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/box/fvelementgeometry.hh:379
friend std::ranges::view auto boundaryFaces(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:443
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/box/fvelementgeometry.hh:375
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/box/fvelementgeometry.hh:564
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:511
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/box/fvelementgeometry.hh:537
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/box/fvelementgeometry.hh:560
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/box/fvelementgeometry.hh:373
void bind(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:503
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element (2^dim for cubes)
Definition discretization/box/fvelementgeometry.hh:385
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/box/fvelementgeometry.hh:232
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/box/fvelementgeometry.hh:169
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/box/fvelementgeometry.hh:293
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:214
const Element & element() const
The bound element.
Definition discretization/box/fvelementgeometry.hh:236
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/box/fvelementgeometry.hh:271
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:110
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/box/fvelementgeometry.hh:263
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvelementgeometry.hh:181
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/box/fvelementgeometry.hh:77
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/box/fvelementgeometry.hh:93
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/box/fvelementgeometry.hh:175
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/box/fvelementgeometry.hh:187
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/box/fvelementgeometry.hh:71
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/box/fvelementgeometry.hh:255
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/box/fvelementgeometry.hh:133
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/box/fvelementgeometry.hh:259
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:123
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/box/fvelementgeometry.hh:305
typename GG::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvelementgeometry.hh:79
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/box/fvelementgeometry.hh:240
void bind(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:206
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvelementgeometry.hh:81
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/box/fvelementgeometry.hh:248
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/box/fvelementgeometry.hh:278
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/box/fvelementgeometry.hh:325
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition discretization/box/fvelementgeometry.hh:244
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/box/fvelementgeometry.hh:336
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/box/fvelementgeometry.hh:151
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/box/fvelementgeometry.hh:315
friend std::ranges::view auto boundaryFaces(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:144
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/box/fvelementgeometry.hh:267
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element (2^dim for cubes)
Definition discretization/box/fvelementgeometry.hh:83
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/box/fvelementgeometry.hh:73
void bindElement(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:223
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/box/fvelementgeometry.hh:99
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/box/fvelementgeometry.hh:89
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:197
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/box/fvelementgeometry.hh:69
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/box/fvelementgeometry.hh:75
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition discretization/box/fvelementgeometry.hh:46
An interpolation point related to a face of an element.
Definition cvfe/interpolationpointdata.hh:109
An interpolation point related to an element that includes global and local positions.
Definition cvfe/interpolationpointdata.hh:31
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition cvfe/interpolationpointdata.hh:82
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition cvfe/interpolationpointdata.hh:60
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
Quadrature rules over sub-control volumes and sub-control volume faces.
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