version 3.8
discretization/facecentered/staggered/fvgridgeometry.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_GRID_GEOMETRY
13#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_GRID_GEOMETRY
14
15#include <memory>
16
17#include <dune/common/rangeutilities.hh>
18#include <dune/grid/common/scsgmapper.hh>
19
23#include <dumux/common/math.hh>
24
29
37
38namespace Dumux {
39
46template<class GridView>
48{
54
55 template<class GridGeometry>
57
58 template<class GridGeometry, bool enableCache>
60
62 {
63 static constexpr auto dim = GridView::Grid::dimension;
64 static constexpr auto numFacesPerElement = dim * 2;
65 static constexpr auto numScvsPerElement = numFacesPerElement;
66 static constexpr auto numLateralScvfsPerScv = 2 * (dim - 1);
68 static constexpr auto minNumScvfsPerElement = numLateralScvfsPerElement // number of lateral faces
69 + numFacesPerElement; // number of central frontal faces
71 + numFacesPerElement; // number of potential frontal faces on boundary
72 };
73};
74
81template<class GridView,
82 bool cachingEnabled = false,
85
92template<class GV, class Traits>
94: public BaseGridGeometry<GV, Traits>
95{
98 using GridIndexType = typename IndexTraits<GV>::GridIndex;
99 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
100 using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex;
101 using Element = typename GV::template Codim<0>::Entity;
102
103 using IntersectionMapper = typename Traits::IntersectionMapper;
104 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
105
106 using Scalar = typename GV::ctype;
107
108 static constexpr auto dim = Traits::StaticInfo::dim;
109 static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement;
110 static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv;
111 static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement;
112 static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement;
113 static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement;
114
115 using ScvfCornerStorage = typename Traits::SubControlVolumeFace::Traits::CornerStorage;
116 using ScvCornerStorage = typename Traits::SubControlVolume::Traits::CornerStorage;
117
118public:
121 static constexpr DiscretizationMethod discMethod{};
122
123 static constexpr bool cachingEnabled = true;
124
128 using LocalView = typename Traits::template LocalView<ThisType, true>;
130 using SubControlVolume = typename Traits::SubControlVolume;
132 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
134 using GridView = GV;
136 using GeometryHelper = typename Traits::GeometryHelper;
138 using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper;
140 using StaticInformation = typename Traits::StaticInfo;
143
145 FaceCenteredStaggeredFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg, const std::string& paramGroup = "")
146 : ParentType(std::move(gg))
147 , intersectionMapper_(this->gridView())
148 {
149 // Check if the overlap size is what we expect
151 DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs at least an overlap of 1 for parallel computations. "
152 << " Set the parameter \"Grid.Overlap\" in the input file.");
153
154 update_();
155 }
156
158 FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "")
159 : FaceCenteredStaggeredFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView), paramGroup)
160 {}
161
163 std::size_t numScv() const
164 { return scvs_.size(); }
165
167 std::size_t numScvf() const
168 { return scvfs_.size(); }
169
171 std::size_t numBoundaryScv() const
172 { return numBoundaryScv_; }
173
175 std::size_t numBoundaryScvf() const
176 { return numBoundaryScvf_; }
177
179 std::size_t numIntersections() const
180 { return intersectionMapper_.numIntersections(); }
181
183 std::size_t numDofs() const
184 { return this->gridView().size(1); }
185
187 void update(const GridView& gridView)
188 {
189 ParentType::update(gridView);
190 update_();
191 }
192
194 void update(GridView&& gridView)
195 {
196 ParentType::update(std::move(gridView));
197 update_();
198 }
199
201 const SubControlVolume& scv(GridIndexType scvIdx) const
202 { return scvs_[scvIdx]; }
203
206 auto scvs(const LocalView& fvGeometry) const
207 {
208 auto begin = scvs_.cbegin() + numScvsPerElement*fvGeometry.elementIndex();
209 const auto end = begin + numScvsPerElement;
210 return Dune::IteratorRange<std::decay_t<decltype(begin)>>(begin, end);
211 }
212
214 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
215 { return scvfs_[scvfIdx]; }
216
218 const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const
219 { return scvfIndicesOfElement_[eIdx]; }
220
225 const ConnectivityMap& connectivityMap() const
226 { return connectivityMap_; }
227
229 bool hasBoundaryScvf(GridIndexType eIdx) const
230 { return hasBoundaryScvf_[eIdx]; }
231
233 const IntersectionMapper& intersectionMapper() const
234 { return intersectionMapper_; }
235
237 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
238 { return periodicFaceMap_.count(dofIdx); }
239
241 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
242 { return periodicFaceMap_.at(dofIdx); }
243
245 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
246 { return periodicFaceMap_; }
247
248private:
249
250 void update_()
251 {
252 // clear containers (necessary after grid refinement)
253 scvs_.clear();
254 scvfs_.clear();
255 scvfIndicesOfElement_.clear();
256 intersectionMapper_.update(this->gridView());
257
258 // determine size of containers
259 const auto numElements = this->gridView().size(0);
260 scvfIndicesOfElement_.resize(numElements);
261 hasBoundaryScvf_.resize(numElements, false);
262
263 outSideBoundaryVolVarIdx_ = 0;
264 numBoundaryScv_ = 0;
265 numBoundaryScvf_ = 0;
266
267 GeometryHelper geometryHelper(this->gridView());
268
269 // get the global scvf indices first
270 GridIndexType numScvfs = 0;
271 for (const auto& element : elements(this->gridView()))
272 {
273 assert(numScvsPerElement == element.subEntities(1));
274
275 for (const auto& intersection : intersections(this->gridView(), element))
276 {
277 // frontal scvf in element center
278 ++numScvfs;
279
280 // lateral scvfs
281 numScvfs += numLateralScvfsPerScv;
282
283 // handle physical domain boundary
284 if (onDomainBoundary_(intersection))
285 {
286 ++numBoundaryScv_; // frontal face
287 numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces
288
289 // frontal scvf at boundary
290 ++numScvfs;
291 }
292 }
293 }
294
295 // allocate memory
296 const auto numScvs = numElements*numScvsPerElement;
297 scvs_.resize(numScvs);
298 scvfs_.reserve(numScvfs);
299
300 // Build the scvs and scv faces
301 std::size_t globalScvfIdx = 0;
302 for (const auto& element : elements(this->gridView()))
303 {
304 const auto eIdx = this->elementMapper().index(element);
305 auto& globalScvfIndices = scvfIndicesOfElement_[eIdx];
306 globalScvfIndices.resize(minNumScvfsPerElement);
307 globalScvfIndices.reserve(maxNumScvfsPerElement);
308
309 auto getGlobalScvIdx = [&](const auto elementIdx, const auto localScvIdx)
310 { return numScvsPerElement*elementIdx + localScvIdx; };
311
312 LocalIntersectionMapper localIsMapper;
313 localIsMapper.update(this->gridView(), element);
314
315 for (const auto& intersection : intersections(this->gridView(), element))
316 {
317 const auto& intersectionUnitOuterNormal = intersection.centerUnitOuterNormal();
318 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
319 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
320
321 const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx);
322 const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
323 const auto localOppositeScvIdx = geometryHelper.localOppositeIdx(localScvIdx);
324 const auto& intersectionGeometry = intersection.geometry();
325 const auto& elementGeometry = element.geometry();
326
327 assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside());
328
329 // handle periodic boundaries
330 if (onPeriodicBoundary_(intersection))
331 {
332 this->setPeriodic();
333
334 const auto& otherElement = intersection.outside();
335
336 SmallLocalIndexType otherIntersectionLocalIdx = 0;
337 bool periodicFaceFound = false;
338
339 for (const auto& otherIntersection : intersections(this->gridView(), otherElement))
340 {
341 if (periodicFaceFound)
342 continue;
343
344 if (Dune::FloatCmp::eq(intersectionUnitOuterNormal*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7))
345 {
346 const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx);
347 periodicFaceMap_[dofIndex] = periodicDofIdx;
348 periodicFaceFound = true;
349 }
350
351 ++otherIntersectionLocalIdx;
352 }
353 }
354
355 // the sub control volume
356 scvs_[globalScvIdx] = SubControlVolume(
357 elementGeometry,
358 intersectionGeometry,
359 globalScvIdx,
360 localScvIdx,
361 dofIndex,
362 Dumux::normalAxis(intersectionUnitOuterNormal),
363 this->elementMapper().index(element),
364 onDomainBoundary_(intersection)
365 );
366
367 // the frontal sub control volume face at the element center
368 scvfs_.emplace_back(elementGeometry,
369 intersectionGeometry,
370 std::array{globalScvIdx, getGlobalScvIdx(eIdx, localOppositeScvIdx)},
371 localScvfIdx,
372 globalScvfIdx,
373 intersectionUnitOuterNormal,
374 SubControlVolumeFace::FaceType::frontal,
375 SubControlVolumeFace::BoundaryType::interior
376 );
377
378 globalScvfIndices[localScvfIdx] = globalScvfIdx++;
379 ++localScvfIdx;
380
381 // the lateral sub control volume faces
382 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx),
383 [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
384 )
385 {
386 const auto& lateralIntersection = geometryHelper.intersection(lateralFacetIndex, element);
387
388 // helper lambda to get the lateral scvf's global inside and outside scv indices
389 const auto globalScvIndicesForLateralFace = [&]
390 {
391 const auto globalOutsideScvIdx = [&]
392 {
393 if (lateralIntersection.neighbor())
394 {
395 const auto parallelElemIdx = this->elementMapper().index(lateralIntersection.outside());
396 return getGlobalScvIdx(parallelElemIdx, localScvIdx);
397 }
398 else if (onDomainBoundary_(lateralIntersection))
399 return numScvs + outSideBoundaryVolVarIdx_++;
400 else
401 return globalScvIdx; // fallback for parallel, won't be used anyway
402 }();
403
404 return std::array{globalScvIdx, globalOutsideScvIdx};
405 }();
406
407 const auto boundaryType = [&]
408 {
409 if (onProcessorBoundary_(lateralIntersection))
410 return SubControlVolumeFace::BoundaryType::processorBoundary;
411 else if (onDomainBoundary_(lateralIntersection))
412 return SubControlVolumeFace::BoundaryType::physicalBoundary;
413 else
414 return SubControlVolumeFace::BoundaryType::interior;
415 }();
416
417 scvfs_.emplace_back(
418 elementGeometry,
419 intersectionGeometry,
420 geometryHelper.facet(lateralFacetIndex, element).geometry(),
421 globalScvIndicesForLateralFace, // TODO higher order
422 localScvfIdx,
423 globalScvfIdx,
424 lateralIntersection.centerUnitOuterNormal(),
425 SubControlVolumeFace::FaceType::lateral,
426 boundaryType
427 );
428
429 globalScvfIndices[localScvfIdx] = globalScvfIdx++;
430 ++localScvfIdx;
431
432 if (onDomainBoundary_(lateralIntersection))
433 {
434 ++numBoundaryScvf_;
435 hasBoundaryScvf_[eIdx] = true;
436 }
437 } // end loop over lateral facets
438
439 } // end first loop over intersections
440
441 // do a second loop over all intersections to add frontal boundary faces
442 int localScvfIdx = minNumScvfsPerElement;
443 for (const auto& intersection : intersections(this->gridView(), element))
444 {
445 // the frontal sub control volume face at a domain boundary (coincides with element face)
446 if (onDomainBoundary_(intersection))
447 {
448 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
449 const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx);
450 ++numBoundaryScvf_;
451
452 // the frontal sub control volume face at the boundary
453 scvfs_.emplace_back(
454 element.geometry(),
455 intersection.geometry(),
456 std::array{globalScvIdx, globalScvIdx}, // TODO outside boundary, periodic, parallel?
457 localScvfIdx,
458 globalScvfIdx,
459 intersection.centerUnitOuterNormal(),
460 SubControlVolumeFace::FaceType::frontal,
461 SubControlVolumeFace::BoundaryType::physicalBoundary
462 );
463
464 globalScvfIndices.push_back(globalScvfIdx);
465 ++globalScvfIdx;
466 ++localScvfIdx;
467 hasBoundaryScvf_[eIdx] = true;
468 }
469 }
470 }
471
472 connectivityMap_.update(*this);
473 }
474
475 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
476 {
477 return !intersection.neighbor() && intersection.boundary();
478 }
479
480 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
481 {
482 return !intersection.neighbor() && !intersection.boundary();
483 }
484
485 bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const
486 {
487 return intersection.boundary() && intersection.neighbor();
488 }
489
490 // mappers
491 ConnectivityMap connectivityMap_;
492 IntersectionMapper intersectionMapper_;
493
494 std::vector<SubControlVolume> scvs_;
495 std::vector<SubControlVolumeFace> scvfs_;
496 GridIndexType numBoundaryScv_;
497 GridIndexType numBoundaryScvf_;
498 GridIndexType outSideBoundaryVolVarIdx_;
499 std::vector<bool> hasBoundaryScvf_;
500
501 std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_;
502
503 // a map for periodic boundary vertices
504 std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_;
505};
506
513template<class GV, class Traits>
515: public BaseGridGeometry<GV, Traits>
516{
519 using GridIndexType = typename IndexTraits<GV>::GridIndex;
520 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
521 using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex;
522 using Element = typename GV::template Codim<0>::Entity;
523
524 using IntersectionMapper = typename Traits::IntersectionMapper;
525 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
526
527 static constexpr auto dim = Traits::StaticInfo::dim;
528 static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement;
529 static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv;
530 static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement;
531 static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement;
532 static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement;
533
534public:
537 static constexpr DiscretizationMethod discMethod{};
538
539 static constexpr bool cachingEnabled = false;
540
544 using LocalView = typename Traits::template LocalView<ThisType, false>;
546 using SubControlVolume = typename Traits::SubControlVolume;
548 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
550 using GridView = GV;
552 using GeometryHelper = typename Traits::GeometryHelper;
554 using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper;
556 using StaticInformation = typename Traits::StaticInfo;
559
561 FaceCenteredStaggeredFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg, const std::string& paramGroup = "")
562 : ParentType(std::move(gg))
563 , intersectionMapper_(this->gridView())
564 {
565 // Check if the overlap size is what we expect
567 DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs at least an overlap of 1 for parallel computations. "
568 << " Set the parameter \"Grid.Overlap\" in the input file.");
569
570 update_();
571 }
572
574 FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "")
575 : FaceCenteredStaggeredFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView), paramGroup)
576 {}
577
579 std::size_t numScv() const
580 { return numScvs_; }
581
583 std::size_t numScvf() const
584 { return numScvf_; }
585
587 std::size_t numBoundaryScv() const
588 { return numBoundaryScv_; }
589
591 std::size_t numBoundaryScvf() const
592 { return numBoundaryScvf_; }
593
595 std::size_t numIntersections() const
596 { return intersectionMapper_.numIntersections(); }
597
599 std::size_t numDofs() const
600 { return this->gridView().size(1); }
601
606 const ConnectivityMap& connectivityMap() const
607 { return connectivityMap_; }
608
610 bool hasBoundaryScvf(GridIndexType eIdx) const
611 { return hasBoundaryScvf_[eIdx]; }
612
614 const IntersectionMapper& intersectionMapper() const
615 { return intersectionMapper_; }
616
618 const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const
619 { return scvfIndicesOfElement_[eIdx]; }
620
622 GridIndexType outsideVolVarIndex(GridIndexType scvfIdx) const
623 { return outsideVolVarIndices_.at(scvfIdx); }
624
626 void update(const GridView& gridView)
627 {
628 ParentType::update(gridView);
629 update_();
630 }
631
633 void update(GridView&& gridView)
634 {
635 ParentType::update(std::move(gridView));
636 update_();
637 }
638
640 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
641 { return periodicFaceMap_.count(dofIdx); }
642
644 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
645 { return periodicFaceMap_.at(dofIdx); }
646
648 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
649 { return periodicFaceMap_; }
650
651private:
652
653 void update_()
654 {
655 intersectionMapper_.update(this->gridView());
656
657 // clear local data
658 numScvf_ = 0;
659 numBoundaryScv_ = 0;
660 numBoundaryScvf_ = 0;
661 hasBoundaryScvf_.clear();
662 scvfIndicesOfElement_.clear();
663 outsideVolVarIndices_.clear();
664
665 // determine size of containers
666 const auto numElements = this->gridView().size(0);
667 scvfIndicesOfElement_.resize(numElements);
668 hasBoundaryScvf_.resize(numElements, false);
669 numScvs_ = numElements*numScvsPerElement;
670
671 GeometryHelper geometryHelper(this->gridView());
672
673 // get the global scv indices first
674 GridIndexType scvfIdx = 0;
675
676 GridIndexType neighborVolVarIdx = numScvs_;
677
678 for (const auto& element : elements(this->gridView()))
679 {
680 const auto eIdx = this->elementMapper().index(element);
681 assert(numScvsPerElement == element.subEntities(1));
682
683 // the element-wise index sets for finite volume geometry
684 auto& globalScvfIndices = scvfIndicesOfElement_[eIdx];
685 globalScvfIndices.reserve(maxNumScvfsPerElement);
686 globalScvfIndices.resize(minNumScvfsPerElement);
687
688 // keep track of frontal boundary scvfs
689 std::size_t numFrontalBoundaryScvfs = 0;
690
691 using LocalIntersectionIndexMapper = FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView>;
692 LocalIntersectionIndexMapper localIsMapper;
693 localIsMapper.update(this->gridView(), element);
694
695 for (const auto& intersection : intersections(this->gridView(), element))
696 {
697 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
698 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
699
700 assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside());
701 // the frontal sub control volume face at the element center
702 globalScvfIndices[localScvfIdx] = scvfIdx++;
703 ++localScvfIdx;
704
705 if constexpr(dim > 1)
706 {
707 // the lateral sub control volume faces
708 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx),
709 [&](auto idx) { return localIsMapper.refToRealIdx(idx) ;})
710 )
711 {
712 if (onDomainBoundary_(geometryHelper.intersection(lateralFacetIndex, element)))
713 {
714 outsideVolVarIndices_[scvfIdx] = neighborVolVarIdx++;
715 ++numBoundaryScvf_;
716 hasBoundaryScvf_[eIdx] = true;
717 }
718
719 globalScvfIndices[localScvfIdx] = scvfIdx++;
720 ++localScvfIdx;
721 }
722 }
723
724 // handle physical domain boundary
725 if (onDomainBoundary_(intersection))
726 {
727 ++numBoundaryScv_; // frontal face
728 numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces
729 ++numFrontalBoundaryScvfs;
730 ++numBoundaryScvf_;
731 hasBoundaryScvf_[eIdx] = true;
732 }
733
734 // handle periodic boundaries
735 if (onPeriodicBoundary_(intersection))
736 {
737 this->setPeriodic();
738
739 const auto& otherElement = intersection.outside();
740
741 SmallLocalIndexType otherIntersectionLocalIdx = 0;
742 bool periodicFaceFound = false;
743
744 for (const auto& otherIntersection : intersections(this->gridView(), otherElement))
745 {
746 if (periodicFaceFound)
747 continue;
748
749 if (Dune::FloatCmp::eq(intersection.centerUnitOuterNormal()*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7))
750 {
751 const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx);
752 const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, localScvIdx);
753 periodicFaceMap_[dofIndex] = periodicDofIdx;
754 periodicFaceFound = true;
755 }
756
757 ++otherIntersectionLocalIdx;
758 }
759 }
760 }
761
762 // add global indices of frontal boundary scvfs last
763 for (std::size_t i = 0; i < numFrontalBoundaryScvfs; ++i)
764 globalScvfIndices.push_back(scvfIdx++);
765 }
766
767 // set number of subcontrolvolume faces
768 numScvf_ = scvfIdx;
769
770 connectivityMap_.update(*this);
771 }
772
773 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
774 {
775 return !intersection.neighbor() && intersection.boundary();
776 }
777
778 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
779 {
780 return !intersection.neighbor() && !intersection.boundary();
781 }
782
783 bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const
784 {
785 return intersection.boundary() && intersection.neighbor();
786 }
787
788 // mappers
789 ConnectivityMap connectivityMap_;
790 IntersectionMapper intersectionMapper_;
791
793 std::size_t numScvs_;
794 std::size_t numScvf_;
795 std::size_t numBoundaryScv_;
796 std::size_t numBoundaryScvf_;
797 std::vector<bool> hasBoundaryScvf_;
798
799 std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_;
800
801 // a map for periodic boundary vertices
802 std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_;
803 std::unordered_map<GridIndexType, GridIndexType> outsideVolVarIndices_;
804};
805
806} // end namespace Dumux
807
808#endif
Base class for grid geometries.
Check the overlap size for different discretization methods.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
defines a standard intersection mapper for mapping of global DOFs assigned to faces....
Definition: intersectionmapper.hh:29
Definition: localintersectionindexmapper.hh:29
Stores the dof indices corresponding to the neighboring scvs that contribute to the derivative calcul...
Definition: facecentered/staggered/connectivitymap.hh:30
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:98
Base class for the finite volume geometry vector for face-centered staggered models This builds up th...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:516
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:548
typename Traits::GeometryHelper GeometryHelper
export the geometry helper type
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:552
FaceCenteredStaggeredFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg, const std::string &paramGroup="")
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:561
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:546
std::size_t numDofs() const
the total number of dofs
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:599
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the d.o.f. on the other side of the periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:644
const ConnectivityMap & connectivityMap() const
Returns the connectivity map of which dofs have derivatives with respect to a given dof.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:606
std::size_t numBoundaryScv() const
The total number of boundary sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:587
const IntersectionMapper & intersectionMapper() const
Return a reference to the intersection mapper.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:614
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:633
FaceCenteredStaggeredFVGridGeometry(const GridView &gridView, const std::string &paramGroup="")
Constructor from gridView.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:574
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a d.o.f. is on a periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:640
typename Traits::StaticInfo StaticInformation
export static information
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:556
const std::vector< GridIndexType > & scvfIndicesOfElement(GridIndexType eIdx) const
Get the global sub control volume face indices of an element.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:618
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:610
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries // TODO rename to periodic dof map in fvassem...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:648
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:583
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:542
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:591
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:579
typename Traits::LocalIntersectionMapper LocalIntersectionMapper
export the local intersection mapper
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:554
GridIndexType outsideVolVarIndex(GridIndexType scvfIdx) const
Get the global sub control volume face indices of an element.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:622
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:626
std::size_t numIntersections() const
The total number of intersections.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:595
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:558
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:544
Base class for the finite volume geometry vector for staggered models This builds up the sub control ...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:95
const ConnectivityMap & connectivityMap() const
Returns the connectivity map of which dofs have derivatives with respect to a given dof.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:225
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:167
const IntersectionMapper & intersectionMapper() const
Return a reference to the intersection mapper.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:233
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:229
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:187
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:132
typename Traits::LocalIntersectionMapper LocalIntersectionMapper
export the local intersection mapper
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:138
typename Traits::GeometryHelper GeometryHelper
export the geometry helper type
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:136
const std::vector< GridIndexType > & scvfIndicesOfElement(GridIndexType eIdx) const
Get the global sub control volume face indices of an element.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:218
FaceCenteredStaggeredFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg, const std::string &paramGroup="")
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:145
typename Traits::StaticInfo StaticInformation
export static information
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:140
std::size_t numDofs() const
the total number of dofs
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:183
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries // TODO rename to periodic dof map in fvassem...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:245
std::size_t numBoundaryScv() const
The total number of boundary sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:171
FaceCenteredStaggeredFVGridGeometry(const GridView &gridView, const std::string &paramGroup="")
Constructor from gridView.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:158
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the d.o.f. on the other side of the periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:241
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:201
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:130
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:128
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:194
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:126
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:142
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:175
std::size_t numIntersections() const
The total number of intersections.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:179
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:214
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a d.o.f. is on a periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:237
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:163
auto scvs(const LocalView &fvGeometry) const
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:206
Base class for the finite volume geometry vector for face-centered staggered models This builds up th...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:84
Face centered staggered geometry helper.
Definition: discretization/facecentered/staggered/geometryhelper.hh:28
Face centered staggered sub control volume face.
Definition: discretization/facecentered/staggered/subcontrolvolumeface.hh:54
Face centered staggered sub control volume.
Definition: discretization/facecentered/staggered/subcontrolvolume.hh:52
Defines the default element and vertex mapper types.
Geometry helper for face-centered staggered scheme.
Helper classes to compute the integration elements.
Dune::Std::detected_or_t< Dumux::BasicGridGeometry< GV, typename T::ElementMapper, typename T::VertexMapper >, Detail::SpecifiesBaseGridGeometry, T > BasicGridGeometry_t
Type of the basic grid geometry implementation used as backend.
Definition: basegridgeometry.hh:42
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
Defines the index types used for grid and local indices.
defines intersection mappers.
Define some often used mathematical functions.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:28
Definition: defaultmappertraits.hh:23
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:62
static constexpr auto numLateralScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:67
static constexpr auto maxNumScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:70
static constexpr auto numScvsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:65
static constexpr auto numFacesPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:64
static constexpr auto dim
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:63
static constexpr auto numLateralScvfsPerScv
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:66
static constexpr auto minNumScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:68
The default traits for the face-center staggered finite volume grid geometry Defines the scv and scvf...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:48
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26