version 3.8
discretization/facecentered/staggered/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//
12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
14
15#include <utility>
16#include <bitset>
17
18#include <dune/common/rangeutilities.hh>
19#include <dune/common/reservedvector.hh>
20
22#include <dune/common/iteratorrange.hh>
27
28namespace Dumux {
29
30#ifndef DOXYGEN
31namespace Detail::FCStaggered {
32
33template<class FVElementGeometry, class SubControlVolume>
34typename SubControlVolume::Traits::Geometry scvGeometry(const FVElementGeometry& fvGeometry,
35 const SubControlVolume& scv)
36{
37 typename SubControlVolume::Traits::CornerStorage corners{};
38 // select the containing element
39 const auto elementGeometry = (scv.elementIndex() != fvGeometry.elementIndex()) ?
40 fvGeometry.element().geometry() :
41 fvGeometry.gridGeometry().element(scv.elementIndex()).geometry();
42
43 const auto center = elementGeometry.center();
44 const auto dofAxis = scv.dofAxis();
45 for (int i = 0; i < corners.size(); ++i)
46 {
47 auto& corner = corners[i];
48
49 // copy the corner of the corresponding element
50 corner = elementGeometry.corner(i);
51
52 // shift the corner such that the scv covers half of the element
53 // (keep the corner positions at the face with the staggered dof)
54 if ((corner[dofAxis] - center[dofAxis]) * scv.directionSign() < 0.0)
55 corner[dofAxis] = center[dofAxis];
56 }
57
58 return {corners.front(), corners.back()};
59}
60
61template<class FVElementGeometry, class SubControlVolumeFace>
62typename SubControlVolumeFace::Traits::Geometry scvfGeometry(const FVElementGeometry& fvGeometry,
63 const SubControlVolumeFace& scvf)
64{
65 const auto normalAxis = scvf.normalAxis();
66 const auto center = scvf.center();
67 const auto shift = scvf.ipGlobal() - center;
68 const auto dofAxis = scvf.isLateral() ? Dumux::normalAxis(shift) : normalAxis;
69 const auto insideElementIndex = (fvGeometry.scv(scvf.insideScvIdx())).elementIndex();
70 const auto elementGeometry = (insideElementIndex != fvGeometry.elementIndex()) ?
71 fvGeometry.element().geometry() :
72 fvGeometry.gridGeometry().element(insideElementIndex).geometry();
73
74 auto corners = std::array{
75 elementGeometry.corner(0),
76 elementGeometry.corner(elementGeometry.corners() - 1)
77 };
78
79 // shift corners to scvf plane and halve lateral faces
80 for (int i = 0; i < corners.size(); ++i)
81 {
82 auto& corner = corners[i];
83 corner[normalAxis] = center[normalAxis];
84 if (scvf.isLateral() && (corner - center)*shift < 0.0)
85 corner[dofAxis] = elementGeometry.center()[dofAxis];
86 }
87
88 auto inPlaneAxes = std::move(std::bitset<SubControlVolumeFace::Traits::dimWorld>{}.set());
89 inPlaneAxes.set(normalAxis, false);
90
91 return {corners[0], corners[1], inPlaneAxes};
92}
93
94} // end namespace Detail::FCStaggered
95#endif // DOXYGEN
96
97template<class GG, bool cachingEnabled>
99
105template<class GG>
106class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true>
107{
108 using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true>;
109 using GridView = typename GG::GridView;
110 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
111 static constexpr auto numScvsPerElement = GG::StaticInformation::numScvsPerElement;
112
113public:
115 using SubControlVolume = typename GG::SubControlVolume;
116 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
117 using Element = typename GridView::template Codim<0>::Entity;
118 using GridGeometry = GG;
119
121 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
122
124 : gridGeometry_(&gridGeometry)
125 {}
126
128 const SubControlVolume& scv(GridIndexType scvIdx) const
129 { return gridGeometry().scv(scvIdx); }
130
132 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
133 { return gridGeometry().scvf(scvfIdx); }
134
137 {
138 assert(scvf.isLateral());
139 const auto otherGlobalIdx = scvfIndices_()[GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex())];
140 return gridGeometry().scvf(otherGlobalIdx);
141
142 }
143
146 {
147 assert(scv.boundary());
148
149 // frontal boundary faces are always stored after the lateral faces
150 auto scvfIter = scvfs(*this, scv).begin();
151 const auto end = scvfs(*this, scv).end();
152 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
153 ++scvfIter;
154
155 assert(scvfIter->isFrontal());
156 assert(scvfIter->boundary());
157 return *scvfIter;
158 }
159
165 friend inline auto
167 { return fvGeometry.gridGeometry().scvs(fvGeometry); }
168
174 friend inline auto
176 {
177 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
179 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(fvGeometry.scvfIndices_().begin(), fvGeometry),
180 ScvfIterator(fvGeometry.scvfIndices_().end(), fvGeometry));
181 }
182
188 friend inline auto
190 {
191 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
193 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
194 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
195 return Dune::IteratorRange<ScvfIterator>(begin, end);
196 }
197
199 std::size_t numScv() const
200 { return numScvsPerElement; }
201
203 std::size_t numScvf() const
204 { return scvfIndices_().size(); }
205
207 bool hasBoundaryScvf() const
208 { return gridGeometry().hasBoundaryScvf(eIdx_); }
209
216 {
217 this->bind(element);
218 return std::move(*this);
219 }
220
222 void bind(const Element& element) &
223 { this->bindElement(element); }
224
231 {
232 this->bindElement(element);
233 return std::move(*this);
234 }
235
237 void bindElement(const Element& element) &
238 {
239 elementPtr_ = &element;
240 eIdx_ = gridGeometry().elementMapper().index(element);
241 }
242
245 {
246 assert(gridGeometry_);
247 return *gridGeometry_;
248 }
249
250 std::size_t elementIndex() const
251 { return eIdx_; }
252
254 const Element& element() const
255 { return *elementPtr_; }
256
259 { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); }
260
263 {
264 const auto& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
265 assert(!lateralOrthogonalScvf.boundary());
266
267 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
268
269 auto outsideFVGeometry = localView(gridGeometry());
270 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
271 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
272
273 for (const auto& otherScvf : scvfs(outsideFVGeometry))
274 if (otherScvf.localIndex() == otherLocalIdx)
275 return otherScvf;
276
277 DUNE_THROW(Dune::InvalidStateException, "No outside scvf found");
278 }
279
281 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
282 {
283 return Detail::FCStaggered::scvGeometry(*this, scv);
284 }
285
287 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
288 {
289 return Detail::FCStaggered::scvfGeometry(*this, scvf);
290 }
291
292private:
293
294 const auto& scvfIndices_() const
295 {
296 return gridGeometry().scvfIndicesOfElement(eIdx_);
297 }
298
299 const Element* elementPtr_;
300 GridIndexType eIdx_;
301 const GridGeometry* gridGeometry_;
302};
303
309template<class GG>
310class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>
311{
312 using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>;
313
314 using GridView = typename GG::GridView;
315 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
316
317 //TODO include assert that checks for quad geometry
318
319 using StaticInfo = typename GG::StaticInformation;
320 static constexpr auto dim = StaticInfo::dim;
321 static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement;
322 static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv;
323 static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement;
324 static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement;
325 static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement;
326
327 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
328
329public:
331 using SubControlVolume = typename GG::SubControlVolume;
332 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
333 using Element = typename GridView::template Codim<0>::Entity;
334 using GridGeometry = GG;
335
337 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
338
340 : gridGeometry_(&gridGeometry)
341 , geometryHelper_(gridGeometry.gridView())
342 {}
343
345 const SubControlVolumeFace& scvf(const GridIndexType scvfIdx) const
346 { return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; }
347
353 friend inline auto
355 {
356 using Iter = typename decltype(fvGeometry.scvfs_)::const_iterator;
357 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
358 }
359
365 friend inline auto
367 {
368 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
370 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
371 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
372 return Dune::IteratorRange<ScvfIterator>(begin, end);
373 }
374
376 const SubControlVolume& scv(const GridIndexType scvIdx) const
377 {
378 if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back())
379 return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)];
380 else
381 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
382 }
383
389 friend inline auto
391 {
392 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
393 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
394 }
395
398 {
399 assert(scvf.isLateral());
400 const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex());
401 return scvfs_[otherLocalIdx];
402 }
403
406 {
407 assert(scv.boundary());
408
409 // frontal boundary faces are always stored after the lateral faces
410 auto scvfIter = scvfs(*this, scv).begin();
411 const auto end = scvfs(*this, scv).end();
412 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
413 ++scvfIter;
414
415 assert(scvfIter->isFrontal());
416 assert(scvfIter->boundary());
417 return *scvfIter;
418 }
419
421 bool hasBoundaryScvf() const
422 { return gridGeometry().hasBoundaryScvf(eIdx_); }
423
425 std::size_t numScv() const
426 { return numScvsPerElement; }
427
429 std::size_t numScvf() const
430 { return scvfIndices_().size(); }
431
438 {
439 this->bind_(element);
440 return std::move(*this);
441 }
442
443 void bind(const Element& element) &
444 { this->bind_(element); }
445
452 {
453 typename GG::LocalIntersectionMapper localIsMapper;
454 localIsMapper.update(gridGeometry().gridView(), element);
455 this->bindElement_(element, localIsMapper);
456 return std::move(*this);
457 }
458
459 void bindElement(const Element& element) &
460 {
461 typename GG::LocalIntersectionMapper localIsMapper;
462 localIsMapper.update(gridGeometry().gridView(), element);
463 this->bindElement_(element, localIsMapper);
464 }
465
466 GridIndexType elementIndex() const
467 { return eIdx_; }
468
470 const Element& element() const
471 { return *elementPtr_; }
472
475 {
476 assert(gridGeometry_);
477 return *gridGeometry_;
478 }
479
482 { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); }
483
487 {
488 const SubControlVolumeFace& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
489 assert(!lateralOrthogonalScvf.boundary());
490
491 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
492
493 auto outsideFVGeometry = localView(gridGeometry());
494 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
495 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
496
497 for (const auto& otherScvf : scvfs(outsideFVGeometry))
498 if (otherScvf.localIndex() == otherLocalIdx)
499 return otherScvf;
500
501 DUNE_THROW(Dune::InvalidStateException, "No outside scvf found");
502 }
503
505 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
506 {
507 return Detail::FCStaggered::scvGeometry(*this, scv);
508 }
509
511 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
512 {
513 return Detail::FCStaggered::scvfGeometry(*this, scvf);
514 }
515
516private:
519 void bind_(const Element& element)
520 {
521 typename GG::LocalIntersectionMapper localIsMapper;
522 localIsMapper.update(gridGeometry().gridView(), element);
523
524 bindElement_(element, localIsMapper);
525 neighborScvIndices_.clear();
526 neighborScvs_.clear();
527
528 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
529 {
530 if (intersection.neighbor())
531 {
532 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
533 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
534 const auto& neighborElement = intersection.outside();
535 const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement);
536 const auto& neighborElementGeometry = neighborElement.geometry();
537
538 // todo: could be done easier?
539 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
540 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement);
541
542 typename GG::LocalIntersectionMapper localNeighborIsMapper;
543 localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement);
544
545 for (const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement))
546 {
547 const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside());
548 if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx)
549 {
550
551 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside());
552 neighborScvs_.push_back(SubControlVolume(
553 neighborElementGeometry,
554 neighborIntersection.geometry(),
555 globalScvIndicesOfNeighborElement[localNeighborScvIdx],
556 localNeighborScvIdx,
557 dofIndex,
558 Dumux::normalAxis(neighborIntersection.centerUnitOuterNormal()),
559 neighborElementIdx,
560 onDomainBoundary_(neighborIntersection)
561 ));
562
563 neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]);
564 }
565 }
566 }
567 }
568 }
569
571 void bindElement_(const Element& element, const typename GG::LocalIntersectionMapper& localIsMapper)
572 {
573 elementPtr_ = &element;
574 eIdx_ = gridGeometry().elementMapper().index(element);
575 std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement);
576 scvfs_.clear();
577 scvfs_.resize(minNumScvfsPerElement);
578
579 const auto& elementGeometry = element.geometry();
580
581 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
582 {
583 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
584 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
585 const auto& intersectionGeometry = intersection.geometry();
586 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
587
588 scvs_[localScvIdx] = SubControlVolume(
589 elementGeometry,
590 intersectionGeometry,
591 scvIndicesOfElement_[localScvIdx],
592 localScvIdx,
593 dofIndex,
594 Dumux::normalAxis(intersection.centerUnitOuterNormal()),
595 eIdx_,
596 onDomainBoundary_(intersection)
597 );
598
599 // the frontal sub control volume face at the element center
600 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
601 scvfs_[localScvfIdx] = SubControlVolumeFace(
602 elementGeometry,
603 intersectionGeometry,
604 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]},
605 localScvfIdx,
606 scvfIndices_()[localScvfIdx],
607 intersection.centerUnitOuterNormal(),
608 SubControlVolumeFace::FaceType::frontal,
609 SubControlVolumeFace::BoundaryType::interior
610 );
611 ++localScvfIdx;
612
613 // the lateral sub control volume faces
614 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx),
615 [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
616 )
617 {
618 const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element);
619 const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element);
620 const auto& lateralFacetGeometry = lateralFacet.geometry();
621
622 // helper lambda to get the lateral scvf's global inside and outside scv indices
623 const auto globalScvIndicesForLateralFace = [&]
624 {
625 const auto globalOutsideScvIdx = [&]
626 {
627 if (lateralIntersection.neighbor())
628 {
629 const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside());
630 // todo: could be done easier?
631 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
632 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement);
633 return globalScvIndicesOfNeighborElement[localScvIdx];
634 }
635 else
636 return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]);
637 }();
638
639 return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx};
640 }();
641
642 const auto boundaryType = [&]
643 {
644 if (onProcessorBoundary_(lateralIntersection))
645 return SubControlVolumeFace::BoundaryType::processorBoundary;
646 else if (onDomainBoundary_(lateralIntersection))
647 return SubControlVolumeFace::BoundaryType::physicalBoundary;
648 else
649 return SubControlVolumeFace::BoundaryType::interior;
650 }();
651
652 scvfs_[localScvfIdx] = SubControlVolumeFace(
653 elementGeometry,
654 intersectionGeometry,
655 lateralFacetGeometry,
656 globalScvIndicesForLateralFace, // TODO higher order
657 localScvfIdx,
658 scvfIndices_()[localScvfIdx],
659 lateralIntersection.centerUnitOuterNormal(),
660 SubControlVolumeFace::FaceType::lateral,
661 boundaryType
662 );
663 ++localScvfIdx;
664 }
665 }
666
667 // do a second loop over all intersections to add frontal boundary faces
668 auto localScvfIdx = minNumScvfsPerElement;
669 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
670 {
671 // the frontal sub control volume face at a domain boundary (coincides with element face)
672 if (onDomainBoundary_(intersection))
673 {
674 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
675 // the frontal sub control volume face at the boundary
676 scvfs_.push_back(SubControlVolumeFace(
677 element.geometry(),
678 intersection.geometry(),
679 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]},
680 localScvfIdx,
681 scvfIndices_()[localScvfIdx],
682 intersection.centerUnitOuterNormal(),
683 SubControlVolumeFace::FaceType::frontal,
684 SubControlVolumeFace::BoundaryType::physicalBoundary
685 ));
686 ++localScvfIdx;
687 }
688 }
689
690 if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{})
691 {
692 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
693 if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement)
694 {
695 // make sure frontal boundary scvfs are sorted correctly
696 std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(),
697 [](const auto& scvfLeft, const auto& scvfRight)
698 { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); }
699 );
700 }
701 }
702 }
703
704 const auto& scvfIndices_() const
705 { return gridGeometry().scvfIndicesOfElement(eIdx_); }
706
707 template<class Entry, class Container>
708 const LocalIndexType findLocalIndex_(const Entry& entry,
709 const Container& container) const
710 {
711 auto it = std::find(container.begin(), container.end(), entry);
712 assert(it != container.end() && "Could not find the entry! Make sure to properly bind this class!");
713 return std::distance(container.begin(), it);
714 }
715
716 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
717 {
718 return !intersection.neighbor() && intersection.boundary();
719 }
720
721 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
722 {
723 return !intersection.neighbor() && !intersection.boundary();
724 }
725
726 Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_;
727
728 Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_;
729 Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_;
730
731 std::array<SubControlVolume, numScvsPerElement> scvs_;
732 std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_;
733
734 const GridGeometry* gridGeometry_;
735 const Element* elementPtr_;
736 GridIndexType eIdx_;
737 typename GridGeometry::GeometryHelper geometryHelper_;
738};
739
740} // end namespace Dumux
741
742#endif
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:311
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:429
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:421
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:354
void bindElement(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:459
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:331
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:437
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:481
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:451
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:334
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:511
GridIndexType elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:466
SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:486
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:333
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:397
const SubControlVolume & scv(const GridIndexType scvIdx) const
Get a half sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:376
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:474
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &g)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:390
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:505
void bind(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:443
const SubControlVolumeFace & scvf(const GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:345
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:470
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:405
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:366
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:339
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:425
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:332
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:107
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:115
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:215
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:230
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:199
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:254
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:244
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:175
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:132
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:287
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:207
void bind(const Element &element) &
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:222
const SubControlVolumeFace & outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Returns the the scvf of neighbor element with the same integration point and unit outer normal.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:262
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:166
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:128
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:116
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:123
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:118
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:258
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:281
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:117
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:189
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:136
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:203
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:145
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:237
std::size_t elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:250
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:98
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:70
Iterators over sub control volume faces of an fv geometry and a given sub control volume.
Definition: scvandscvfiterators.hh:110
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
static std::size_t normalAxis(const Vector &v)
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: normalaxis.hh:26
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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