version 3.7
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 (std::size_t 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 scvfs_.emplace_back(MpfaHelper(),
341 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
342 is,
343 vIdxGlobal,
344 vIdxLocal,
345 scvfIdx,
346 eIdx,
347 outsideScvIndices,
348 q,
349 boundary);
350
351 // insert the scvf data into the dual grid index set
352 dualIdSet[vIdxGlobal].insert(scvfs_.back());
353
354 // increment scvf counter
355 scvfIdx++;
356 }
357
358 // for network grids, clear outside indices to not make a second scvf on that facet
359 if (dim < dimWorld)
360 outsideIndices[indexInInside].clear();
361 }
362
363 // create sub control volume for this element
364 scvs_[eIdx] = SubControlVolume(std::move(elementGeometry), eIdx);
365
366 // Save the scvf indices belonging to this scv to build up fv element geometries fast
367 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
368 }
369
370 // Make the flip index set for network and surface grids
371 flipScvfIndices_.resize(scvfs_.size());
372 for (const auto& scvf : scvfs_)
373 {
374 if (scvf.boundary())
375 continue;
376
377 const auto numOutsideScvs = scvf.numOutsideScvs();
378 const auto vIdxGlobal = scvf.vertexIndex();
379 const auto insideScvIdx = scvf.insideScvIdx();
380
381 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
382 for (std::size_t i = 0; i < numOutsideScvs; ++i)
383 {
384 const auto outsideScvIdx = scvf.outsideScvIdx(i);
385 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
386 {
387 const auto& outsideScvf = this->scvf(outsideScvfIndex);
388 if (outsideScvf.vertexIndex() == vIdxGlobal &&
389 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
390 {
391 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
392 // there is always only one flip face in an outside element
393 break;
394 }
395 }
396 }
397 }
398
399 // building the geometries has finished
400 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
401
402 // Initialize the grid interaction volume index sets
403 timer.reset();
404 ivIndexSets_.update(*this, std::move(dualIdSet));
405 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
406
407 // build the connectivity map for an efficient assembly
408 timer.reset();
409 connectivityMap_.update(*this);
410 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
411 }
412
413 // connectivity map for efficient assembly
414 ConnectivityMap connectivityMap_;
415
416 // the finite volume grid geometries
417 std::vector<SubControlVolume> scvs_;
418 std::vector<SubControlVolumeFace> scvfs_;
419
420 // containers storing the global data
421 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
422 std::vector<bool> secondaryInteractionVolumeVertices_;
423 GridIndexType numBoundaryScvf_;
424 std::vector<bool> hasBoundaryScvf_;
425
426 // needed for embedded surface and network grids (dim < dimWorld)
427 FlipScvfIndexSet flipScvfIndices_;
428
429 // The grid interaction volume index set
430 GridIVIndexSets ivIndexSets_;
431
432 // Indicator function on where to use the secondary IVs
433 SecondaryIvIndicatorType secondaryIvIndicator_;
434};
435
443template<class GV, class Traits>
444class CCMpfaFVGridGeometry<GV, Traits, false>
445: public BaseGridGeometry<GV, Traits>
446{
449
450 static constexpr int dim = GV::dimension;
451 static constexpr int dimWorld = GV::dimensionworld;
452
453 using Element = typename GV::template Codim<0>::Entity;
454 using Vertex = typename GV::template Codim<dim>::Entity;
455 using Intersection = typename GV::Intersection;
456 using GridIndexType = typename IndexTraits<GV>::GridIndex;
457 using CoordScalar = typename GV::ctype;
458
459 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
460
461public:
463 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
465 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
467 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
468
470 using LocalView = typename Traits::template LocalView<ThisType, false>;
472 using SubControlVolume = typename Traits::SubControlVolume;
474 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
478 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
480 using DofMapper = typename Traits::ElementMapper;
482 using GridView = GV;
484 using MpfaHelper = typename Traits::template MpfaHelper<ThisType>;
485
488 static constexpr DiscretizationMethod discMethod{};
489
491 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
492
494 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
495
499 : ParentType(gridView)
500 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
501 { return is.boundary() || isBranching; } )
502 {
503 checkOverlapSizeCCMpfa(gridView);
504 update_();
505 }
506
508 CCMpfaFVGridGeometry(const GridView& gridView, const SecondaryIvIndicatorType& indicator)
509 : ParentType(gridView)
510 , secondaryIvIndicator_(indicator)
511 {
512 checkOverlapSizeCCMpfa(gridView);
513 update_();
514 }
515
518 const DofMapper& dofMapper() const
519 { return this->elementMapper(); }
520
522 std::size_t numScv() const
523 { return numScvs_; }
524
526 std::size_t numScvf() const
527 { return numScvf_; }
528
530 std::size_t numBoundaryScvf() const
531 { return numBoundaryScvf_; }
532
534 std::size_t numDofs() const
535 { return this->gridView().size(0); }
536
539 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary, bool> = 0>
540 bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
541 { return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
542
545 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary, bool> = 0>
546 constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
547 { return false; }
548
550 bool isGhostVertex(const Vertex& v) const
551 { return isGhostVertex_[this->vertexMapper().index(v)]; }
552
554 bool isGhostVertex(GridIndexType vIdxGlobal) const
555 { return isGhostVertex_[vIdxGlobal]; }
556
557
559 void update(const GridView& gridView)
560 {
561 ParentType::update(gridView);
562 update_();
563 }
564
566 void update(GridView&& gridView)
567 {
568 ParentType::update(std::move(gridView));
569 update_();
570 }
571
574 { return MpfaHelper(); }
575
577 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
578 { return scvfIndicesOfScv_[scvIdx]; }
579
581 const std::vector<ScvfOutsideGridIndexStorage>& neighborVolVarIndices(GridIndexType scvIdx) const
582 { return neighborVolVarIndices_[scvIdx]; }
583
586 const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
587 { return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
588
591 { return flipScvfIndices_; }
592
596 { return connectivityMap_; }
597
600 { return ivIndexSets_; }
601
602private:
603
604 void update_()
605 {
606
607 // stop the time required for the update
608 Dune::Timer timer;
609
610 // resize containers
611 numScvs_ = numDofs();
612 scvfIndicesOfScv_.resize(numScvs_);
613 neighborVolVarIndices_.resize(numScvs_);
614
615 // Some methods require to use a second type of interaction volume, e.g.
616 // around vertices on the boundary or branching points (surface grids)
617 const auto numVert = this->gridView().size(dim);
618 secondaryInteractionVolumeVertices_.assign(numVert, false);
619
620 // find vertices on processor boundaries HERE!!
621 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
622
623 // instantiate the dual grid index set (to be used for construction of interaction volumes)
624 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
625
626 // keep track of boundary scvfs and scvf vertex indices in order to set up flip scvf index set
627 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
628 std::vector<bool> scvfIsOnBoundary;
629 std::vector<GridIndexType> scvfVertexIndex;
630 scvfIsOnBoundary.reserve(maxNumScvfs);
631 scvfVertexIndex.reserve(maxNumScvfs);
632
633 // Build the SCVs and SCV faces
634 numScvf_ = 0;
635 numBoundaryScvf_ = 0;
636 for (const auto& element : elements(this->gridView()))
637 {
638 const auto eIdx = this->elementMapper().index(element);
639
640 // the element geometry
641 auto elementGeometry = element.geometry();
642
643 // the element-wise index sets for finite volume geometry
644 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
645 std::vector<GridIndexType> scvfsIndexSet;
646 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
647 scvfsIndexSet.reserve(numLocalFaces);
648 neighborVolVarIndexSet.reserve(numLocalFaces);
649
650 // for network grids there might be multiple intersections with the same geometryInInside
651 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
652 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
653 if (dim < dimWorld)
654 {
655 outsideIndices.resize(element.subEntities(1));
656 for (const auto& intersection : intersections(this->gridView(), element))
657 {
658 if (intersection.neighbor())
659 {
660 auto nIdx = this->elementMapper().index( intersection.outside() );
661 outsideIndices[intersection.indexInInside()].push_back(nIdx);
662 }
663 }
664 }
665
666 // construct the sub control volume faces
667 for (const auto& is : intersections(this->gridView(), element))
668 {
669 const auto indexInInside = is.indexInInside();
670 const bool boundary = is.boundary();
671 const bool neighbor = is.neighbor();
672
673 // for surface grids, skip the rest if handled already
674 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
675 continue;
676
677 // if outside level > inside level, use the outside element in the following
678 const bool useNeighbor = neighbor && is.outside().level() > element.level();
679 const auto& e = useNeighbor ? is.outside() : element;
680 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
681 const auto eg = e.geometry();
682 const auto refElement = referenceElement(eg);
683
684 // evaluate if vertices on this intersection use primary/secondary IVs
685 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
686 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
687
688 // make the scv faces belonging to each corner of the intersection
689 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
690 {
691 // get the global vertex index the scv face is connected to
692 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
693 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
694
695 // do not build scvfs connected to a processor boundary
696 if (isGhostVertex_[vIdxGlobal])
697 continue;
698
699 // if this vertex is tagged to use the secondary IVs, store info
700 if (usesSecondaryIV)
701 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
702
703 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
704 const auto& outsideScvIndices = [&] ()
705 {
706 if (!boundary)
707 return dim == dimWorld ?
708 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
709 outsideIndices[indexInInside];
710 else
711 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
712 } ();
713
714 // insert the scvf data into the dual grid index set
715 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
716
717 // store information on the scv face
718 scvfsIndexSet.push_back(numScvf_++);
719 scvfIsOnBoundary.push_back(boundary);
720 scvfVertexIndex.push_back(vIdxGlobal);
721 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
722 }
723
724 // for network grids, clear outside indices to not make a second scvf on that facet
725 if (dim < dimWorld)
726 outsideIndices[indexInInside].clear();
727 }
728
729 // store the sets of indices in the data container
730 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
731 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
732 }
733
734 // Make the flip scvf index set
735 flipScvfIndices_.resize(numScvf_);
736 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
737 {
738 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
739 for (unsigned int i = 0; i < scvfIndices.size(); ++i)
740 {
741 // boundary scvf have no flip scvfs
742 if (scvfIsOnBoundary[ scvfIndices[i] ])
743 continue;
744
745 const auto scvfIdx = scvfIndices[i];
746 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
747 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
748
749 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
750 for (unsigned int j = 0; j < numOutsideScvs; ++j)
751 {
752 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
753 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
754 for (unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
755 {
756 const auto outsideScvfIndex = outsideScvfIndices[k];
757 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
758 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
759 if (outsideScvfVertexIndex == vIdxGlobal &&
760 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
761 {
762 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
763 // there is always only one flip face in an outside element
764 break;
765 }
766 }
767 }
768 }
769 }
770
771 // building the geometries has finished
772 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
773
774 // Initialize the grid interaction volume index sets
775 timer.reset();
776 ivIndexSets_.update(*this, std::move(dualIdSet));
777 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
778
779 // build the connectivity map for an efficient assembly
780 timer.reset();
781 connectivityMap_.update(*this);
782 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
783 }
784
785 // connectivity map for efficient assembly
786 ConnectivityMap connectivityMap_;
787
788 // containers storing the global data
789 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
790 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
791 std::vector<bool> secondaryInteractionVolumeVertices_;
792 std::vector<bool> isGhostVertex_;
793 GridIndexType numScvs_;
794 GridIndexType numScvf_;
795 GridIndexType numBoundaryScvf_;
796
797 // needed for embedded surface and network grids (dim < dimWorld)
798 FlipScvfIndexSet flipScvfIndices_;
799
800 // The grid interaction volume index set
801 GridIVIndexSets ivIndexSets_;
802
803 // Indicator function on where to use the secondary IVs
804 SecondaryIvIndicatorType secondaryIvIndicator_;
805};
806
807} // end namespace Dumux
808
809#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:446
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:467
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:474
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:478
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:554
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:586
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:566
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:599
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:476
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:465
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:581
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:508
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:590
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:498
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:522
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:577
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:463
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:472
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:530
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:480
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:573
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:518
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:595
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:550
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:484
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:470
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:540
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:546
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:526
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:559
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:534
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