3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_GRID_GEOMETRY
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_GRID_GEOMETRY
26
27#include <memory>
28
29#include <dune/common/rangeutilities.hh>
30#include <dune/grid/common/scsgmapper.hh>
31
35#include <dumux/common/math.hh>
36
41
49
50namespace Dumux {
51
58template<class GridView>
60{
61 using ElementMapper = Dune::SingleCodimSingleGeomTypeMapper<GridView, 0/*codim*/>;
67
68 template<class GridGeometry>
70
71 template<class GridGeometry, bool enableCache>
73
75 {
76 static constexpr auto dim = GridView::Grid::dimension;
77 static constexpr auto numFacesPerElement = dim * 2;
78 static constexpr auto numScvsPerElement = numFacesPerElement;
79 static constexpr auto numLateralScvfsPerScv = 2 * (dim - 1);
81 static constexpr auto minNumScvfsPerElement = numLateralScvfsPerElement // number of lateral faces
82 + numFacesPerElement; // number of central frontal faces
84 + numFacesPerElement; // number of potential frontal faces on boundary
85 };
86};
87
94template<class GridView,
95 bool cachingEnabled = false,
98
105template<class GV, class Traits>
107: public BaseGridGeometry<GV, Traits>
108{
111 using GridIndexType = typename IndexTraits<GV>::GridIndex;
112 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
113 using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex;
114 using Element = typename GV::template Codim<0>::Entity;
115
116 using IntersectionMapper = typename Traits::IntersectionMapper;
117 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
118
119 using Scalar = typename GV::ctype;
120
121 static constexpr auto dim = Traits::StaticInfo::dim;
122 static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement;
123 static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv;
124 static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement;
125 static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement;
126 static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement;
127
128 using ScvfCornerStorage = typename Traits::SubControlVolumeFace::Traits::CornerStorage;
129 using ScvCornerStorage = typename Traits::SubControlVolume::Traits::CornerStorage;
130
131public:
134 static constexpr DiscretizationMethod discMethod{};
135
136 static constexpr bool cachingEnabled = true;
137
139 using LocalView = typename Traits::template LocalView<ThisType, true>;
141 using SubControlVolume = typename Traits::SubControlVolume;
143 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
145 using GridView = GV;
147 using GeometryHelper = typename Traits::GeometryHelper;
149 using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper;
151 using StaticInformation = typename Traits::StaticInfo;
154
156 FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "")
157 : ParentType(gridView)
158 , intersectionMapper_(gridView)
159 {
160 // Check if the overlap size is what we expect
162 DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs at least an overlap of 1 for parallel computations. "
163 << " Set the parameter \"Grid.Overlap\" in the input file.");
164
165 update_();
166 }
167
169 std::size_t numScv() const
170 { return scvs_.size(); }
171
173 std::size_t numScvf() const
174 { return scvfs_.size(); }
175
177 std::size_t numBoundaryScv() const
178 { return numBoundaryScv_; }
179
181 std::size_t numBoundaryScvf() const
182 { return numBoundaryScvf_; }
183
185 std::size_t numIntersections() const
186 { return intersectionMapper_.numIntersections(); }
187
189 std::size_t numDofs() const
190 { return this->gridView().size(1); }
191
193 void update(const GridView& gridView)
194 {
195 ParentType::update(gridView);
196 update_();
197 }
198
200 void update(GridView&& gridView)
201 {
202 ParentType::update(std::move(gridView));
203 update_();
204 }
205
207 const SubControlVolume& scv(GridIndexType scvIdx) const
208 { return scvs_[scvIdx]; }
209
212 auto scvs(const LocalView& fvGeometry) const
213 {
214 auto begin = scvs_.cbegin() + numScvsPerElement*fvGeometry.elementIndex();
215 const auto end = begin + numScvsPerElement;
216 return Dune::IteratorRange<std::decay_t<decltype(begin)>>(begin, end);
217 }
218
220 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
221 { return scvfs_[scvfIdx]; }
222
224 const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const
225 { return scvfIndicesOfElement_[eIdx]; }
226
231 const ConnectivityMap& connectivityMap() const
232 { return connectivityMap_; }
233
235 bool hasBoundaryScvf(GridIndexType eIdx) const
236 { return hasBoundaryScvf_[eIdx]; }
237
239 const IntersectionMapper& intersectionMapper() const
240 { return intersectionMapper_; }
241
243 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
244 { return periodicFaceMap_.count(dofIdx); }
245
247 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
248 { return periodicFaceMap_.at(dofIdx); }
249
251 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
252 { return periodicFaceMap_; }
253
254private:
255
256 void update_()
257 {
258 // clear containers (necessary after grid refinement)
259 scvs_.clear();
260 scvfs_.clear();
261 scvfIndicesOfElement_.clear();
262 intersectionMapper_.update(this->gridView());
263
264 // determine size of containers
265 const auto numElements = this->gridView().size(0);
266 scvfIndicesOfElement_.resize(numElements);
267 hasBoundaryScvf_.resize(numElements, false);
268
269 outSideBoundaryVolVarIdx_ = 0;
270 numBoundaryScv_ = 0;
271 numBoundaryScvf_ = 0;
272
273 GeometryHelper geometryHelper(this->gridView());
274
275 // get the global scvf indices first
276 GridIndexType numScvfs = 0;
277 for (const auto& element : elements(this->gridView()))
278 {
279 assert(numScvsPerElement == element.subEntities(1));
280
281 for (const auto& intersection : intersections(this->gridView(), element))
282 {
283 // frontal scvf in element center
284 ++numScvfs;
285
286 // lateral scvfs
287 numScvfs += numLateralScvfsPerScv;
288
289 // handle physical domain boundary
290 if (onDomainBoundary_(intersection))
291 {
292 ++numBoundaryScv_; // frontal face
293 numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces
294
295 // frontal scvf at boundary
296 ++numScvfs;
297 }
298 }
299 }
300
301 // allocate memory
302 const auto numScvs = numElements*numScvsPerElement;
303 scvs_.resize(numScvs);
304 scvfs_.reserve(numScvfs);
305
306 // Build the scvs and scv faces
307 std::size_t globalScvfIdx = 0;
308 for (const auto& element : elements(this->gridView()))
309 {
310 const auto eIdx = this->elementMapper().index(element);
311 auto& globalScvfIndices = scvfIndicesOfElement_[eIdx];
312 globalScvfIndices.resize(minNumScvfsPerElement);
313 globalScvfIndices.reserve(maxNumScvfsPerElement);
314
315 auto getGlobalScvIdx = [&](const auto elementIdx, const auto localScvIdx)
316 { return numScvsPerElement*elementIdx + localScvIdx; };
317
318 LocalIntersectionMapper localIsMapper;
319 localIsMapper.update(this->gridView(), element);
320
321 for (const auto& intersection : intersections(this->gridView(), element))
322 {
323 const auto& intersectionUnitOuterNormal = intersection.centerUnitOuterNormal();
324 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
325 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
326
327 const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx);
328 const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
329 const auto localOppositeScvIdx = geometryHelper.localOppositeIdx(localScvIdx);
330 const auto& intersectionGeometry = intersection.geometry();
331 const auto& elementGeometry = element.geometry();
332
333 assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside());
334
335 // handle periodic boundaries
336 if (onPeriodicBoundary_(intersection))
337 {
338 this->setPeriodic();
339
340 const auto& otherElement = intersection.outside();
341
342 SmallLocalIndexType otherIntersectionLocalIdx = 0;
343 bool periodicFaceFound = false;
344
345 for (const auto& otherIntersection : intersections(this->gridView(), otherElement))
346 {
347 if (periodicFaceFound)
348 continue;
349
350 if (Dune::FloatCmp::eq(intersectionUnitOuterNormal*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7))
351 {
352 const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx);
353 periodicFaceMap_[dofIndex] = periodicDofIdx;
354 periodicFaceFound = true;
355 }
356
357 ++otherIntersectionLocalIdx;
358 }
359 }
360
361 // the sub control volume
362 scvs_[globalScvIdx] = SubControlVolume(
363 elementGeometry,
364 intersectionGeometry,
365 globalScvIdx,
366 localScvIdx,
367 dofIndex,
368 Dumux::normalAxis(intersectionUnitOuterNormal),
369 this->elementMapper().index(element),
370 onDomainBoundary_(intersection)
371 );
372
373 // the frontal sub control volume face at the element center
374 scvfs_.emplace_back(elementGeometry,
375 intersectionGeometry,
376 std::array{globalScvIdx, getGlobalScvIdx(eIdx, localOppositeScvIdx)},
377 localScvfIdx,
378 globalScvfIdx,
379 intersectionUnitOuterNormal,
380 SubControlVolumeFace::FaceType::frontal,
381 SubControlVolumeFace::BoundaryType::interior
382 );
383
384 globalScvfIndices[localScvfIdx] = globalScvfIdx++;
385 ++localScvfIdx;
386
387 // the lateral sub control volume faces
388 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx),
389 [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
390 )
391 {
392 const auto& lateralIntersection = geometryHelper.intersection(lateralFacetIndex, element);
393
394 // helper lambda to get the lateral scvf's global inside and outside scv indices
395 const auto globalScvIndicesForLateralFace = [&]
396 {
397 const auto globalOutsideScvIdx = [&]
398 {
399 if (lateralIntersection.neighbor())
400 {
401 const auto parallelElemIdx = this->elementMapper().index(lateralIntersection.outside());
402 return getGlobalScvIdx(parallelElemIdx, localScvIdx);
403 }
404 else if (onDomainBoundary_(lateralIntersection))
405 return numScvs + outSideBoundaryVolVarIdx_++;
406 else
407 return globalScvIdx; // fallback for parallel, won't be used anyway
408 }();
409
410 return std::array{globalScvIdx, globalOutsideScvIdx};
411 }();
412
413 const auto boundaryType = [&]
414 {
415 if (onProcessorBoundary_(lateralIntersection))
416 return SubControlVolumeFace::BoundaryType::processorBoundary;
417 else if (onDomainBoundary_(lateralIntersection))
418 return SubControlVolumeFace::BoundaryType::physicalBoundary;
419 else
420 return SubControlVolumeFace::BoundaryType::interior;
421 }();
422
423 scvfs_.emplace_back(
424 elementGeometry,
425 intersectionGeometry,
426 geometryHelper.facet(lateralFacetIndex, element).geometry(),
427 globalScvIndicesForLateralFace, // TODO higher order
428 localScvfIdx,
429 globalScvfIdx,
430 lateralIntersection.centerUnitOuterNormal(),
431 SubControlVolumeFace::FaceType::lateral,
432 boundaryType
433 );
434
435 globalScvfIndices[localScvfIdx] = globalScvfIdx++;
436 ++localScvfIdx;
437
438 if (onDomainBoundary_(lateralIntersection))
439 {
440 ++numBoundaryScvf_;
441 hasBoundaryScvf_[eIdx] = true;
442 }
443 } // end loop over lateral facets
444
445 } // end first loop over intersections
446
447 // do a second loop over all intersections to add frontal boundary faces
448 int localScvfIdx = minNumScvfsPerElement;
449 for (const auto& intersection : intersections(this->gridView(), element))
450 {
451 // the frontal sub control volume face at a domain boundary (coincides with element face)
452 if (onDomainBoundary_(intersection))
453 {
454 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
455 const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx);
456 ++numBoundaryScvf_;
457
458 // the frontal sub control volume face at the boundary
459 scvfs_.emplace_back(
460 element.geometry(),
461 intersection.geometry(),
462 std::array{globalScvIdx, globalScvIdx}, // TODO outside boundary, periodic, parallel?
463 localScvfIdx,
464 globalScvfIdx,
465 intersection.centerUnitOuterNormal(),
466 SubControlVolumeFace::FaceType::frontal,
467 SubControlVolumeFace::BoundaryType::physicalBoundary
468 );
469
470 globalScvfIndices.push_back(globalScvfIdx);
471 ++globalScvfIdx;
472 ++localScvfIdx;
473 hasBoundaryScvf_[eIdx] = true;
474 }
475 }
476 }
477
478 connectivityMap_.update(*this);
479 }
480
481 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
482 {
483 return !intersection.neighbor() && intersection.boundary();
484 }
485
486 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
487 {
488 return !intersection.neighbor() && !intersection.boundary();
489 }
490
491 bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const
492 {
493 return intersection.boundary() && intersection.neighbor();
494 }
495
496 // mappers
497 ConnectivityMap connectivityMap_;
498 IntersectionMapper intersectionMapper_;
499
500 std::vector<SubControlVolume> scvs_;
501 std::vector<SubControlVolumeFace> scvfs_;
502 GridIndexType numBoundaryScv_;
503 GridIndexType numBoundaryScvf_;
504 GridIndexType outSideBoundaryVolVarIdx_;
505 std::vector<bool> hasBoundaryScvf_;
506
507 std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_;
508
509 // a map for periodic boundary vertices
510 std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_;
511};
512
519template<class GV, class Traits>
521: public BaseGridGeometry<GV, Traits>
522{
525 using GridIndexType = typename IndexTraits<GV>::GridIndex;
526 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
527 using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex;
528 using Element = typename GV::template Codim<0>::Entity;
529
530 using IntersectionMapper = typename Traits::IntersectionMapper;
531 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
532
533 static constexpr auto dim = Traits::StaticInfo::dim;
534 static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement;
535 static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv;
536 static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement;
537 static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement;
538 static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement;
539
540public:
543 static constexpr DiscretizationMethod discMethod{};
544
545 static constexpr bool cachingEnabled = false;
546
548 using LocalView = typename Traits::template LocalView<ThisType, false>;
550 using SubControlVolume = typename Traits::SubControlVolume;
552 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
554 using GridView = GV;
556 using GeometryHelper = typename Traits::GeometryHelper;
558 using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper;
560 using StaticInformation = typename Traits::StaticInfo;
563
565 FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "")
566 : ParentType(gridView)
567 , intersectionMapper_(gridView)
568 {
569 // Check if the overlap size is what we expect
571 DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs at least an overlap of 1 for parallel computations. "
572 << " Set the parameter \"Grid.Overlap\" in the input file.");
573
574 update_();
575 }
576
578 std::size_t numScv() const
579 { return numScvs_; }
580
582 std::size_t numScvf() const
583 { return numScvf_; }
584
586 std::size_t numBoundaryScv() const
587 { return numBoundaryScv_; }
588
590 std::size_t numBoundaryScvf() const
591 { return numBoundaryScvf_; }
592
594 std::size_t numIntersections() const
595 { return intersectionMapper_.numIntersections(); }
596
598 std::size_t numDofs() const
599 { return this->gridView().size(1); }
600
605 const ConnectivityMap& connectivityMap() const
606 { return connectivityMap_; }
607
609 bool hasBoundaryScvf(GridIndexType eIdx) const
610 { return hasBoundaryScvf_[eIdx]; }
611
613 const IntersectionMapper& intersectionMapper() const
614 { return intersectionMapper_; }
615
617 const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const
618 { return scvfIndicesOfElement_[eIdx]; }
619
621 GridIndexType outsideVolVarIndex(GridIndexType scvfIdx) const
622 { return outsideVolVarIndices_.at(scvfIdx); }
623
625 void update(const GridView& gridView)
626 {
627 ParentType::update(gridView);
628 update_();
629 }
630
632 void update(GridView&& gridView)
633 {
634 ParentType::update(std::move(gridView));
635 update_();
636 }
637
639 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
640 { return periodicFaceMap_.count(dofIdx); }
641
643 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
644 { return periodicFaceMap_.at(dofIdx); }
645
647 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
648 { return periodicFaceMap_; }
649
650private:
651
652 void update_()
653 {
654 intersectionMapper_.update(this->gridView());
655
656 // clear local data
657 numScvf_ = 0;
658 numBoundaryScv_ = 0;
659 numBoundaryScvf_ = 0;
660 hasBoundaryScvf_.clear();
661 scvfIndicesOfElement_.clear();
662 outsideVolVarIndices_.clear();
663
664 // determine size of containers
665 const auto numElements = this->gridView().size(0);
666 scvfIndicesOfElement_.resize(numElements);
667 hasBoundaryScvf_.resize(numElements, false);
668 numScvs_ = numElements*numScvsPerElement;
669
670 GeometryHelper geometryHelper(this->gridView());
671
672 // get the global scv indices first
673 GridIndexType scvfIdx = 0;
674
675 GridIndexType neighborVolVarIdx = numScvs_;
676
677 for (const auto& element : elements(this->gridView()))
678 {
679 const auto eIdx = this->elementMapper().index(element);
680 assert(numScvsPerElement == element.subEntities(1));
681
682 // the element-wise index sets for finite volume geometry
683 auto& globalScvfIndices = scvfIndicesOfElement_[eIdx];
684 globalScvfIndices.reserve(maxNumScvfsPerElement);
685 globalScvfIndices.resize(minNumScvfsPerElement);
686
687 // keep track of frontal boundary scvfs
688 std::size_t numFrontalBoundaryScvfs = 0;
689
690 using LocalIntersectionIndexMapper = FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView>;
691 LocalIntersectionIndexMapper localIsMapper;
692 localIsMapper.update(this->gridView(), element);
693
694 for (const auto& intersection : intersections(this->gridView(), element))
695 {
696 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
697 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
698
699 assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside());
700 // the frontal sub control volume face at the element center
701 globalScvfIndices[localScvfIdx] = scvfIdx++;
702 ++localScvfIdx;
703
704 if constexpr(dim > 1)
705 {
706 // the lateral sub control volume faces
707 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx),
708 [&](auto idx) { return localIsMapper.refToRealIdx(idx) ;})
709 )
710 {
711 if (onDomainBoundary_(geometryHelper.intersection(lateralFacetIndex, element)))
712 {
713 outsideVolVarIndices_[scvfIdx] = neighborVolVarIdx++;
714 ++numBoundaryScvf_;
715 hasBoundaryScvf_[eIdx] = true;
716 }
717
718 globalScvfIndices[localScvfIdx] = scvfIdx++;
719 ++localScvfIdx;
720 }
721 }
722
723 // handle physical domain boundary
724 if (onDomainBoundary_(intersection))
725 {
726 ++numBoundaryScv_; // frontal face
727 numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces
728 ++numFrontalBoundaryScvfs;
729 ++numBoundaryScvf_;
730 hasBoundaryScvf_[eIdx] = true;
731 }
732
733 // handle periodic boundaries
734 if (onPeriodicBoundary_(intersection))
735 {
736 this->setPeriodic();
737
738 const auto& otherElement = intersection.outside();
739
740 SmallLocalIndexType otherIntersectionLocalIdx = 0;
741 bool periodicFaceFound = false;
742
743 for (const auto& otherIntersection : intersections(this->gridView(), otherElement))
744 {
745 if (periodicFaceFound)
746 continue;
747
748 if (Dune::FloatCmp::eq(intersection.centerUnitOuterNormal()*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7))
749 {
750 const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx);
751 const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, localScvIdx);
752 periodicFaceMap_[dofIndex] = periodicDofIdx;
753 periodicFaceFound = true;
754 }
755
756 ++otherIntersectionLocalIdx;
757 }
758 }
759 }
760
761 // add global indices of frontal boundary scvfs last
762 for (std::size_t i = 0; i < numFrontalBoundaryScvfs; ++i)
763 globalScvfIndices.push_back(scvfIdx++);
764 }
765
766 // set number of subcontrolvolume faces
767 numScvf_ = scvfIdx;
768
769 connectivityMap_.update(*this);
770 }
771
772 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
773 {
774 return !intersection.neighbor() && intersection.boundary();
775 }
776
777 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
778 {
779 return !intersection.neighbor() && !intersection.boundary();
780 }
781
782 bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const
783 {
784 return intersection.boundary() && intersection.neighbor();
785 }
786
787 // mappers
788 ConnectivityMap connectivityMap_;
789 IntersectionMapper intersectionMapper_;
790
792 std::size_t numScvs_;
793 std::size_t numScvf_;
794 std::size_t numBoundaryScv_;
795 std::size_t numBoundaryScvf_;
796 std::vector<bool> hasBoundaryScvf_;
797
798 std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_;
799
800 // a map for periodic boundary vertices
801 std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_;
802 std::unordered_map<GridIndexType, GridIndexType> outsideVolVarIndices_;
803};
804
805} // end namespace Dumux
806
807#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
defines intersection mappers.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
Base class for grid geometries.
Check the overlap size for different discretization methods.
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
static std::size_t normalAxis(const Vector &v)
Definition: normalaxis.hh:39
Definition: defaultmappertraits.hh:35
EM ElementMapper
Definition: defaultmappertraits.hh:36
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
defines a standard intersection mapper for mapping of global DOFs assigned to faces....
Definition: intersectionmapper.hh:41
Base class for all finite volume grid geometries.
Definition: basegridgeometry.hh:51
GV GridView
export the grid view type
Definition: basegridgeometry.hh:66
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:40
Stores the dof indices corresponding to the neighboring scvs that contribute to the derivative calcul...
Definition: facecentered/staggered/connectivitymap.hh:42
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:42
The default traits for the face-center staggered finite volume grid geometry Defines the scv and scvf...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:60
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:75
static constexpr auto numLateralScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:80
static constexpr auto maxNumScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:83
static constexpr auto numScvsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:78
static constexpr auto numFacesPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:77
static constexpr auto dim
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:76
static constexpr auto numLateralScvfsPerScv
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:79
static constexpr auto minNumScvfsPerElement
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:81
Base class for the finite volume geometry vector for face-centered staggered models This builds up th...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:97
Base class for the finite volume geometry vector for staggered models This builds up the sub control ...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:108
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:231
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:173
const IntersectionMapper & intersectionMapper() const
Return a reference to the intersection mapper.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:239
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:235
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:193
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:143
typename Traits::LocalIntersectionMapper LocalIntersectionMapper
export the local intersection mapper
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:149
typename Traits::GeometryHelper GeometryHelper
export the geometry helper type
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:147
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:224
typename Traits::StaticInfo StaticInformation
export static information
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:151
std::size_t numDofs() const
the total number of dofs
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:189
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:251
std::size_t numBoundaryScv() const
The total number of boundary sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:177
FaceCenteredStaggeredFVGridGeometry(const GridView &gridView, const std::string &paramGroup="")
Constructor.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:156
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:247
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:207
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:141
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:139
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:200
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:153
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:181
std::size_t numIntersections() const
The total number of intersections.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:185
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:220
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a d.o.f. is on a periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:243
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:169
auto scvs(const LocalView &fvGeometry) const
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:212
Base class for the finite volume geometry vector for face-centered staggered models This builds up th...
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:522
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:552
typename Traits::GeometryHelper GeometryHelper
export the geometry helper type
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:556
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:550
std::size_t numDofs() const
the total number of dofs
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:598
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:643
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:605
std::size_t numBoundaryScv() const
The total number of boundary sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:586
const IntersectionMapper & intersectionMapper() const
Return a reference to the intersection mapper.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:613
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:632
FaceCenteredStaggeredFVGridGeometry(const GridView &gridView, const std::string &paramGroup="")
Constructor.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:565
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a d.o.f. is on a periodic boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:639
typename Traits::StaticInfo StaticInformation
export static information
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:560
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:617
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:609
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:647
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:582
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:590
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:578
typename Traits::LocalIntersectionMapper LocalIntersectionMapper
export the local intersection mapper
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:558
GridIndexType outsideVolVarIndex(GridIndexType scvfIdx) const
Get the global sub control volume face indices of an element.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:621
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:625
std::size_t numIntersections() const
The total number of intersections.
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:594
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/facecentered/staggered/fvgridgeometry.hh:562
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:548
Definition: discretization/facecentered/staggered/geometryhelper.hh:37
Definition: localintersectionindexmapper.hh:41
Face centered staggered sub control volume.
Definition: discretization/facecentered/staggered/subcontrolvolume.hh:64
Face centered staggered sub control volume face.
Definition: discretization/facecentered/staggered/subcontrolvolumeface.hh:67
Geometry helper for face-centered staggered scheme.