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