version 3.10-dev
discretization/cellcentered/mpfa/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//
14#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
16
17#include <utility>
18
25
26namespace Dumux {
27
39template<class GridView, class Traits, bool enableCache>
41
43template<class GridView>
44void checkOverlapSizeCCMpfa(const GridView& gridView)
45{
46 // Check if the overlap size is what we expect
48 DUNE_THROW(Dune::InvalidStateException, "The ccmpfa discretization method needs at least an overlap of 1 for parallel computations. "
49 << " Set the parameter \"Grid.Overlap\" in the input file.");
50}
51
58template<class GV, class Traits>
59class CCMpfaFVGridGeometry<GV, Traits, true>
60: public BaseGridGeometry<GV, Traits>
61{
64
65 static constexpr int dim = GV::dimension;
66 static constexpr int dimWorld = GV::dimensionworld;
67
68 using Element = typename GV::template Codim<0>::Entity;
69 using Vertex = typename GV::template Codim<dim>::Entity;
70 using Intersection = typename GV::Intersection;
71 using GridIndexType = typename IndexTraits<GV>::GridIndex;
72 using CoordScalar = typename GV::ctype;
73
74 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
75
76public:
78 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
80 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
82 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
83
85 using LocalView = typename Traits::template LocalView<ThisType, true>;
87 using SubControlVolume = typename Traits::SubControlVolume;
89 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
93 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
95 using DofMapper = typename Traits::ElementMapper;
97 using GridView = GV;
99 using MpfaHelper = typename Traits::template MpfaHelper<ThisType>;
100
103 static constexpr DiscretizationMethod discMethod{};
104
106 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
107
109 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
110
114 : ParentType(gridView)
115 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
116 { return is.boundary() || isBranching; } )
117 {
118 checkOverlapSizeCCMpfa(gridView);
119 update_();
120 }
121
123 CCMpfaFVGridGeometry(const GridView& gridView, const SecondaryIvIndicatorType& indicator)
124 : ParentType(gridView)
125 , secondaryIvIndicator_(indicator)
126 {
127 checkOverlapSizeCCMpfa(gridView);
128 update_();
129 }
130
133 const DofMapper& dofMapper() const
134 { return this->elementMapper(); }
135
137 std::size_t numScv() const
138 { return scvs_.size(); }
139
141 std::size_t numScvf() const
142 { return scvfs_.size(); }
143
145 std::size_t numBoundaryScvf() const
146 { return numBoundaryScvf_; }
147
149 std::size_t numDofs() const
150 { return this->gridView().size(0); }
151
154 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary, bool> = 0>
155 bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
156 { return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
157
160 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary, bool> = 0>
161 constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
162 { return false; }
163
164
166 void update(const GridView& gridView)
167 {
168 ParentType::update(gridView);
169 update_();
170 }
171
173 void update(GridView&& gridView)
174 {
175 ParentType::update(std::move(gridView));
176 update_();
177 }
178
181 { return MpfaHelper(); }
182
184 const SubControlVolume& scv(GridIndexType scvIdx) const
185 { return scvs_[scvIdx]; }
186
188 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
189 { return scvfs_[scvfIdx]; }
190
194 { return connectivityMap_; }
195
198 { return ivIndexSets_; }
199
201 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
202 { return scvfIndicesOfScv_[scvIdx]; }
203
206 { return flipScvfIndices_; }
207
210 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
211 { return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]]; }
212
214 bool hasBoundaryScvf(GridIndexType eIdx) const
215 { return hasBoundaryScvf_[eIdx]; }
216
217private:
218
219 void update_()
220 {
221 // stop the time required for the update
222 Dune::Timer timer;
223
224 // clear scvfs container
225 scvfs_.clear();
226
227 // determine the number of geometric entities
228 const auto numVert = this->gridView().size(dim);
229 const auto numScvs = numDofs();
230 std::size_t numScvf = MpfaHelper::getGlobalNumScvf(this->gridView());
231
232 // resize containers
233 scvs_.resize(numScvs);
234 scvfs_.reserve(numScvf);
235 scvfIndicesOfScv_.resize(numScvs);
236 hasBoundaryScvf_.assign(numScvs, false);
237
238 // Some methods require to use a second type of interaction volume, e.g.
239 // around vertices on the boundary or branching points (surface grids)
240 secondaryInteractionVolumeVertices_.assign(numVert, false);
241
242 // find vertices on processor boundaries
243 const auto isGhostVertex = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
244
245 // instantiate the dual grid index set (to be used for construction of interaction volumes)
246 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
247
248 // Build the SCVs and SCV faces
249 GridIndexType scvfIdx = 0;
250 numBoundaryScvf_ = 0;
251 for (const auto& element : elements(this->gridView()))
252 {
253 const auto eIdx = this->elementMapper().index(element);
254
255 // the element geometry
256 auto elementGeometry = element.geometry();
257
258 // The local scvf index set
259 std::vector<GridIndexType> scvfIndexSet;
260 scvfIndexSet.reserve(MpfaHelper::getNumLocalScvfs(elementGeometry.type()));
261
262 // for network grids there might be multiple intersection with the same geometryInInside
263 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
264 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
265 if (dim < dimWorld)
266 {
267 outsideIndices.resize(element.subEntities(1));
268 for (const auto& intersection : intersections(this->gridView(), element))
269 {
270 if (intersection.neighbor())
271 {
272 const auto nIdx = this->elementMapper().index( intersection.outside() );
273 outsideIndices[intersection.indexInInside()].push_back(nIdx);
274 }
275 }
276 }
277
278 // construct the sub control volume faces
279 for (const auto& is : intersections(this->gridView(), element))
280 {
281 const auto indexInInside = is.indexInInside();
282 const bool boundary = is.boundary();
283 const bool neighbor = is.neighbor();
284
285 if (boundary)
286 hasBoundaryScvf_[eIdx] = true;
287
288 // for surface grids, skip the rest if handled already
289 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
290 continue;
291
292 // if outside level > inside level, use the outside element in the following
293 const bool useNeighbor = neighbor && is.outside().level() > element.level();
294 const auto& e = useNeighbor ? is.outside() : element;
295 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
296 const auto eg = e.geometry();
297 const auto refElement = referenceElement(eg);
298
299 // Set up a container with all relevant positions for scvf corner computation
300 const auto numCorners = is.geometry().corners();
301 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
302 refElement,
303 indexInElement,
304 numCorners);
305
306 // evaluate if vertices on this intersection use primary/secondary IVs
307 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
308 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
309
310 // make the scv faces belonging to each corner of the intersection
311 for (int c = 0; c < numCorners; ++c)
312 {
313 // get the global vertex index the scv face is connected to
314 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
315 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
316
317 // do not build scvfs connected to a processor boundary
318 if (isGhostVertex[vIdxGlobal])
319 continue;
320
321 // if this vertex is tagged to use the secondary IVs, store info
322 if (usesSecondaryIV)
323 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
324
325 // the quadrature point parameterizarion to be used on scvfs
326 static const auto q = getParam<CoordScalar>("MPFA.Q");
327
328 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
329 const auto& outsideScvIndices = [&] ()
330 {
331 if (!boundary)
332 return dim == dimWorld ?
333 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
334 outsideIndices[indexInInside];
335 else
336 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs) + numBoundaryScvf_++});
337 } ();
338
339 scvfIndexSet.push_back(scvfIdx);
340 typename SubControlVolumeFace::FacetInfo facetInfo{
341 this->elementMapper().index(e),
342 useNeighbor ? is.indexInOutside() : is.indexInInside(),
343 c
344 };
345 scvfs_.emplace_back(MpfaHelper(),
346 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
347 is,
348 std::move(facetInfo),
349 vIdxGlobal,
350 vIdxLocal,
351 scvfIdx,
352 eIdx,
353 outsideScvIndices,
354 q,
355 boundary);
356
357 // insert the scvf data into the dual grid index set
358 dualIdSet[vIdxGlobal].insert(scvfs_.back());
359
360 // increment scvf counter
361 scvfIdx++;
362 }
363
364 // for network grids, clear outside indices to not make a second scvf on that facet
365 if (dim < dimWorld)
366 outsideIndices[indexInInside].clear();
367 }
368
369 // create sub control volume for this element
370 scvs_[eIdx] = SubControlVolume(std::move(elementGeometry), eIdx);
371
372 // Save the scvf indices belonging to this scv to build up fv element geometries fast
373 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
374 }
375
376 // Make the flip index set for network and surface grids
377 flipScvfIndices_.resize(scvfs_.size());
378 for (const auto& scvf : scvfs_)
379 {
380 if (scvf.boundary())
381 continue;
382
383 const auto numOutsideScvs = scvf.numOutsideScvs();
384 const auto vIdxGlobal = scvf.vertexIndex();
385 const auto insideScvIdx = scvf.insideScvIdx();
386
387 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
388 for (std::size_t i = 0; i < numOutsideScvs; ++i)
389 {
390 const auto outsideScvIdx = scvf.outsideScvIdx(i);
391 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
392 {
393 const auto& outsideScvf = this->scvf(outsideScvfIndex);
394 if (outsideScvf.vertexIndex() == vIdxGlobal &&
395 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
396 {
397 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
398 // there is always only one flip face in an outside element
399 break;
400 }
401 }
402 }
403 }
404
405 // building the geometries has finished
406 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
407
408 // Initialize the grid interaction volume index sets
409 timer.reset();
410 ivIndexSets_.update(*this, std::move(dualIdSet));
411 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
412
413 // build the connectivity map for an efficient assembly
414 timer.reset();
415 connectivityMap_.update(*this);
416 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
417 }
418
419 // connectivity map for efficient assembly
420 ConnectivityMap connectivityMap_;
421
422 // the finite volume grid geometries
423 std::vector<SubControlVolume> scvs_;
424 std::vector<SubControlVolumeFace> scvfs_;
425
426 // containers storing the global data
427 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
428 std::vector<bool> secondaryInteractionVolumeVertices_;
429 GridIndexType numBoundaryScvf_;
430 std::vector<bool> hasBoundaryScvf_;
431
432 // needed for embedded surface and network grids (dim < dimWorld)
433 FlipScvfIndexSet flipScvfIndices_;
434
435 // The grid interaction volume index set
436 GridIVIndexSets ivIndexSets_;
437
438 // Indicator function on where to use the secondary IVs
439 SecondaryIvIndicatorType secondaryIvIndicator_;
440};
441
449template<class GV, class Traits>
450class CCMpfaFVGridGeometry<GV, Traits, false>
451: public BaseGridGeometry<GV, Traits>
452{
455
456 static constexpr int dim = GV::dimension;
457 static constexpr int dimWorld = GV::dimensionworld;
458
459 using Element = typename GV::template Codim<0>::Entity;
460 using Vertex = typename GV::template Codim<dim>::Entity;
461 using Intersection = typename GV::Intersection;
462 using GridIndexType = typename IndexTraits<GV>::GridIndex;
463 using CoordScalar = typename GV::ctype;
464
465 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
466
467public:
469 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
471 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
473 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
474
476 using LocalView = typename Traits::template LocalView<ThisType, false>;
478 using SubControlVolume = typename Traits::SubControlVolume;
480 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
484 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
486 using DofMapper = typename Traits::ElementMapper;
488 using GridView = GV;
490 using MpfaHelper = typename Traits::template MpfaHelper<ThisType>;
491
494 static constexpr DiscretizationMethod discMethod{};
495
497 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
498
500 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
501
505 : ParentType(gridView)
506 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
507 { return is.boundary() || isBranching; } )
508 {
509 checkOverlapSizeCCMpfa(gridView);
510 update_();
511 }
512
514 CCMpfaFVGridGeometry(const GridView& gridView, const SecondaryIvIndicatorType& indicator)
515 : ParentType(gridView)
516 , secondaryIvIndicator_(indicator)
517 {
518 checkOverlapSizeCCMpfa(gridView);
519 update_();
520 }
521
524 const DofMapper& dofMapper() const
525 { return this->elementMapper(); }
526
528 std::size_t numScv() const
529 { return numScvs_; }
530
532 std::size_t numScvf() const
533 { return numScvf_; }
534
536 std::size_t numBoundaryScvf() const
537 { return numBoundaryScvf_; }
538
540 std::size_t numDofs() const
541 { return this->gridView().size(0); }
542
545 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary, bool> = 0>
546 bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
547 { return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
548
551 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary, bool> = 0>
552 constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
553 { return false; }
554
556 bool isGhostVertex(const Vertex& v) const
557 { return isGhostVertex_[this->vertexMapper().index(v)]; }
558
560 bool isGhostVertex(GridIndexType vIdxGlobal) const
561 { return isGhostVertex_[vIdxGlobal]; }
562
563
565 void update(const GridView& gridView)
566 {
567 ParentType::update(gridView);
568 update_();
569 }
570
572 void update(GridView&& gridView)
573 {
574 ParentType::update(std::move(gridView));
575 update_();
576 }
577
580 { return MpfaHelper(); }
581
583 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
584 { return scvfIndicesOfScv_[scvIdx]; }
585
587 const std::vector<ScvfOutsideGridIndexStorage>& neighborVolVarIndices(GridIndexType scvIdx) const
588 { return neighborVolVarIndices_[scvIdx]; }
589
592 const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
593 { return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
594
597 { return flipScvfIndices_; }
598
602 { return connectivityMap_; }
603
606 { return ivIndexSets_; }
607
608private:
609
610 void update_()
611 {
612
613 // stop the time required for the update
614 Dune::Timer timer;
615
616 // resize containers
617 numScvs_ = numDofs();
618 scvfIndicesOfScv_.resize(numScvs_);
619 neighborVolVarIndices_.resize(numScvs_);
620
621 // Some methods require to use a second type of interaction volume, e.g.
622 // around vertices on the boundary or branching points (surface grids)
623 const auto numVert = this->gridView().size(dim);
624 secondaryInteractionVolumeVertices_.assign(numVert, false);
625
626 // find vertices on processor boundaries HERE!!
627 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
628
629 // instantiate the dual grid index set (to be used for construction of interaction volumes)
630 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
631
632 // keep track of boundary scvfs and scvf vertex indices in order to set up flip scvf index set
633 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
634 std::vector<bool> scvfIsOnBoundary;
635 std::vector<GridIndexType> scvfVertexIndex;
636 scvfIsOnBoundary.reserve(maxNumScvfs);
637 scvfVertexIndex.reserve(maxNumScvfs);
638
639 // Build the SCVs and SCV faces
640 numScvf_ = 0;
641 numBoundaryScvf_ = 0;
642 for (const auto& element : elements(this->gridView()))
643 {
644 const auto eIdx = this->elementMapper().index(element);
645
646 // the element geometry
647 auto elementGeometry = element.geometry();
648
649 // the element-wise index sets for finite volume geometry
650 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
651 std::vector<GridIndexType> scvfsIndexSet;
652 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
653 scvfsIndexSet.reserve(numLocalFaces);
654 neighborVolVarIndexSet.reserve(numLocalFaces);
655
656 // for network grids there might be multiple intersections with the same geometryInInside
657 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
658 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
659 if (dim < dimWorld)
660 {
661 outsideIndices.resize(element.subEntities(1));
662 for (const auto& intersection : intersections(this->gridView(), element))
663 {
664 if (intersection.neighbor())
665 {
666 auto nIdx = this->elementMapper().index( intersection.outside() );
667 outsideIndices[intersection.indexInInside()].push_back(nIdx);
668 }
669 }
670 }
671
672 // construct the sub control volume faces
673 for (const auto& is : intersections(this->gridView(), element))
674 {
675 const auto indexInInside = is.indexInInside();
676 const bool boundary = is.boundary();
677 const bool neighbor = is.neighbor();
678
679 // for surface grids, skip the rest if handled already
680 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
681 continue;
682
683 // if outside level > inside level, use the outside element in the following
684 const bool useNeighbor = neighbor && is.outside().level() > element.level();
685 const auto& e = useNeighbor ? is.outside() : element;
686 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
687 const auto eg = e.geometry();
688 const auto refElement = referenceElement(eg);
689
690 // evaluate if vertices on this intersection use primary/secondary IVs
691 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
692 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
693
694 // make the scv faces belonging to each corner of the intersection
695 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
696 {
697 // get the global vertex index the scv face is connected to
698 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
699 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
700
701 // do not build scvfs connected to a processor boundary
702 if (isGhostVertex_[vIdxGlobal])
703 continue;
704
705 // if this vertex is tagged to use the secondary IVs, store info
706 if (usesSecondaryIV)
707 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
708
709 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
710 const auto& outsideScvIndices = [&] ()
711 {
712 if (!boundary)
713 return dim == dimWorld ?
714 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
715 outsideIndices[indexInInside];
716 else
717 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
718 } ();
719
720 // insert the scvf data into the dual grid index set
721 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
722
723 // store information on the scv face
724 scvfsIndexSet.push_back(numScvf_++);
725 scvfIsOnBoundary.push_back(boundary);
726 scvfVertexIndex.push_back(vIdxGlobal);
727 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
728 }
729
730 // for network grids, clear outside indices to not make a second scvf on that facet
731 if (dim < dimWorld)
732 outsideIndices[indexInInside].clear();
733 }
734
735 // store the sets of indices in the data container
736 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
737 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
738 }
739
740 // Make the flip scvf index set
741 flipScvfIndices_.resize(numScvf_);
742 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
743 {
744 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
745 for (unsigned int i = 0; i < scvfIndices.size(); ++i)
746 {
747 // boundary scvf have no flip scvfs
748 if (scvfIsOnBoundary[ scvfIndices[i] ])
749 continue;
750
751 const auto scvfIdx = scvfIndices[i];
752 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
753 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
754
755 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
756 for (unsigned int j = 0; j < numOutsideScvs; ++j)
757 {
758 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
759 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
760 for (unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
761 {
762 const auto outsideScvfIndex = outsideScvfIndices[k];
763 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
764 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
765 if (outsideScvfVertexIndex == vIdxGlobal &&
766 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
767 {
768 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
769 // there is always only one flip face in an outside element
770 break;
771 }
772 }
773 }
774 }
775 }
776
777 // building the geometries has finished
778 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
779
780 // Initialize the grid interaction volume index sets
781 timer.reset();
782 ivIndexSets_.update(*this, std::move(dualIdSet));
783 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
784
785 // build the connectivity map for an efficient assembly
786 timer.reset();
787 connectivityMap_.update(*this);
788 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
789 }
790
791 // connectivity map for efficient assembly
792 ConnectivityMap connectivityMap_;
793
794 // containers storing the global data
795 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
796 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
797 std::vector<bool> secondaryInteractionVolumeVertices_;
798 std::vector<bool> isGhostVertex_;
799 GridIndexType numScvs_;
800 GridIndexType numScvf_;
801 GridIndexType numBoundaryScvf_;
802
803 // needed for embedded surface and network grids (dim < dimWorld)
804 FlipScvfIndexSet flipScvfIndices_;
805
806 // The grid interaction volume index set
807 GridIVIndexSets ivIndexSets_;
808
809 // Indicator function on where to use the secondary IVs
810 SecondaryIvIndicatorType secondaryIvIndicator_;
811};
812
813} // end namespace Dumux
814
815#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
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:452
std::function< bool(const Element &, const Intersection &, bool)> SecondaryIvIndicatorType
export the type to be used for indicators where to use the secondary ivs
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:473
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:480
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:484
bool isGhostVertex(GridIndexType vIdxGlobal) const
Returns true if the vertex (index) lies on a processor boundary inside a ghost element.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:560
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:592
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:572
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:605
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:482
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:471
const std::vector< ScvfOutsideGridIndexStorage > & neighborVolVarIndices(GridIndexType scvIdx) const
Returns the neighboring vol var indices for each scvf contained in an scv.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:587
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:514
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:596
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:504
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:528
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Returns the sub control volume face indices of an scv by global index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:583
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:469
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:478
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:536
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:486
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:579
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:524
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:601
bool isGhostVertex(const Vertex &v) const
Returns true if a given vertex lies on a processor boundary inside a ghost element.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:556
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:490
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:476
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:546
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:552
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:532
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:565
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:540
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:61
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:95
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:89
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:141
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:99
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:210
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Get the sub control volume face indices of an scv by global index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:201
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:91
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:155
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:78
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:193
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:145
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:87
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:80
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:149
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:85
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:166
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:173
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:188
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:137
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume index set class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:197
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:205
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:123
std::function< bool(const Element &, const Intersection &, bool)> SecondaryIvIndicatorType
export the type to be used for indicators where to use the secondary ivs
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:82
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:93
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:133
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:214
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:113
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:180
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:184
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:161
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:40
Helper classes to compute the integration elements.
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
void checkOverlapSizeCCMpfa(const GridView &gridView)
check the overlap size for parallel computations
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:44
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:28
Definition: method.hh:33
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26