3.6-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
176
178 void update(const GridView& gridView)
179 {
180 ParentType::update(gridView);
181 update_();
182 }
183
185 void update(GridView&& gridView)
186 {
187 ParentType::update(std::move(gridView));
188 update_();
189 }
190
193 { return MpfaHelper(); }
194
196 const SubControlVolume& scv(GridIndexType scvIdx) const
197 { return scvs_[scvIdx]; }
198
200 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
201 { return scvfs_[scvfIdx]; }
202
206 { return connectivityMap_; }
207
210 { return ivIndexSets_; }
211
213 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
214 { return scvfIndicesOfScv_[scvIdx]; }
215
218 { return flipScvfIndices_; }
219
222 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
223 { return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]]; }
224
226 bool hasBoundaryScvf(GridIndexType eIdx) const
227 { return hasBoundaryScvf_[eIdx]; }
228
229private:
230
231 void update_()
232 {
233 // stop the time required for the update
234 Dune::Timer timer;
235
236 // clear scvfs container
237 scvfs_.clear();
238
239 // determine the number of geometric entities
240 const auto numVert = this->gridView().size(dim);
241 const auto numScvs = numDofs();
242 std::size_t numScvf = MpfaHelper::getGlobalNumScvf(this->gridView());
243
244 // resize containers
245 scvs_.resize(numScvs);
246 scvfs_.reserve(numScvf);
247 scvfIndicesOfScv_.resize(numScvs);
248 hasBoundaryScvf_.assign(numScvs, false);
249
250 // Some methods require to use a second type of interaction volume, e.g.
251 // around vertices on the boundary or branching points (surface grids)
252 secondaryInteractionVolumeVertices_.assign(numVert, false);
253
254 // find vertices on processor boundaries
255 const auto isGhostVertex = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
256
257 // instantiate the dual grid index set (to be used for construction of interaction volumes)
258 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
259
260 // Build the SCVs and SCV faces
261 GridIndexType scvfIdx = 0;
262 numBoundaryScvf_ = 0;
263 for (const auto& element : elements(this->gridView()))
264 {
265 const auto eIdx = this->elementMapper().index(element);
266
267 // the element geometry
268 auto elementGeometry = element.geometry();
269
270 // The local scvf index set
271 std::vector<GridIndexType> scvfIndexSet;
272 scvfIndexSet.reserve(MpfaHelper::getNumLocalScvfs(elementGeometry.type()));
273
274 // for network grids there might be multiple intersection with the same geometryInInside
275 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
276 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
277 if (dim < dimWorld)
278 {
279 outsideIndices.resize(element.subEntities(1));
280 for (const auto& intersection : intersections(this->gridView(), element))
281 {
282 if (intersection.neighbor())
283 {
284 const auto nIdx = this->elementMapper().index( intersection.outside() );
285 outsideIndices[intersection.indexInInside()].push_back(nIdx);
286 }
287 }
288 }
289
290 // construct the sub control volume faces
291 for (const auto& is : intersections(this->gridView(), element))
292 {
293 const auto indexInInside = is.indexInInside();
294 const bool boundary = is.boundary();
295 const bool neighbor = is.neighbor();
296
297 if (boundary)
298 hasBoundaryScvf_[eIdx] = true;
299
300 // for surface grids, skip the rest if handled already
301 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
302 continue;
303
304 // if outside level > inside level, use the outside element in the following
305 const bool useNeighbor = neighbor && is.outside().level() > element.level();
306 const auto& e = useNeighbor ? is.outside() : element;
307 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
308 const auto eg = e.geometry();
309 const auto refElement = referenceElement(eg);
310
311 // Set up a container with all relevant positions for scvf corner computation
312 const auto numCorners = is.geometry().corners();
313 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
314 refElement,
315 indexInElement,
316 numCorners);
317
318 // evaluate if vertices on this intersection use primary/secondary IVs
319 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
320 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
321
322 // make the scv faces belonging to each corner of the intersection
323 for (std::size_t c = 0; c < numCorners; ++c)
324 {
325 // get the global vertex index the scv face is connected to
326 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
327 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
328
329 // do not build scvfs connected to a processor boundary
330 if (isGhostVertex[vIdxGlobal])
331 continue;
332
333 // if this vertex is tagged to use the secondary IVs, store info
334 if (usesSecondaryIV)
335 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
336
337 // the quadrature point parameterizarion to be used on scvfs
338 static const auto q = getParam<CoordScalar>("MPFA.Q");
339
340 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
341 const auto& outsideScvIndices = [&] ()
342 {
343 if (!boundary)
344 return dim == dimWorld ?
345 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
346 outsideIndices[indexInInside];
347 else
348 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs) + numBoundaryScvf_++});
349 } ();
350
351 scvfIndexSet.push_back(scvfIdx);
352 scvfs_.emplace_back(MpfaHelper(),
353 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
354 is,
355 vIdxGlobal,
356 vIdxLocal,
357 scvfIdx,
358 eIdx,
359 outsideScvIndices,
360 q,
361 boundary);
362
363 // insert the scvf data into the dual grid index set
364 dualIdSet[vIdxGlobal].insert(scvfs_.back());
365
366 // increment scvf counter
367 scvfIdx++;
368 }
369
370 // for network grids, clear outside indices to not make a second scvf on that facet
371 if (dim < dimWorld)
372 outsideIndices[indexInInside].clear();
373 }
374
375 // create sub control volume for this element
376 scvs_[eIdx] = SubControlVolume(std::move(elementGeometry), eIdx);
377
378 // Save the scvf indices belonging to this scv to build up fv element geometries fast
379 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
380 }
381
382 // Make the flip index set for network and surface grids
383 flipScvfIndices_.resize(scvfs_.size());
384 for (const auto& scvf : scvfs_)
385 {
386 if (scvf.boundary())
387 continue;
388
389 const auto numOutsideScvs = scvf.numOutsideScvs();
390 const auto vIdxGlobal = scvf.vertexIndex();
391 const auto insideScvIdx = scvf.insideScvIdx();
392
393 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
394 for (std::size_t i = 0; i < numOutsideScvs; ++i)
395 {
396 const auto outsideScvIdx = scvf.outsideScvIdx(i);
397 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
398 {
399 const auto& outsideScvf = this->scvf(outsideScvfIndex);
400 if (outsideScvf.vertexIndex() == vIdxGlobal &&
401 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
402 {
403 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
404 // there is always only one flip face in an outside element
405 break;
406 }
407 }
408 }
409 }
410
411 // building the geometries has finished
412 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
413
414 // Initialize the grid interaction volume index sets
415 timer.reset();
416 ivIndexSets_.update(*this, std::move(dualIdSet));
417 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
418
419 // build the connectivity map for an efficient assembly
420 timer.reset();
421 connectivityMap_.update(*this);
422 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
423 }
424
425 // connectivity map for efficient assembly
426 ConnectivityMap connectivityMap_;
427
428 // the finite volume grid geometries
429 std::vector<SubControlVolume> scvs_;
430 std::vector<SubControlVolumeFace> scvfs_;
431
432 // containers storing the global data
433 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
434 std::vector<bool> secondaryInteractionVolumeVertices_;
435 GridIndexType numBoundaryScvf_;
436 std::vector<bool> hasBoundaryScvf_;
437
438 // needed for embedded surface and network grids (dim < dimWorld)
439 FlipScvfIndexSet flipScvfIndices_;
440
441 // The grid interaction volume index set
442 GridIVIndexSets ivIndexSets_;
443
444 // Indicator function on where to use the secondary IVs
445 SecondaryIvIndicatorType secondaryIvIndicator_;
446};
447
455template<class GV, class Traits>
456class CCMpfaFVGridGeometry<GV, Traits, false>
457: public BaseGridGeometry<GV, Traits>
458{
461
462 static constexpr int dim = GV::dimension;
463 static constexpr int dimWorld = GV::dimensionworld;
464
465 using Element = typename GV::template Codim<0>::Entity;
466 using Vertex = typename GV::template Codim<dim>::Entity;
467 using Intersection = typename GV::Intersection;
468 using GridIndexType = typename IndexTraits<GV>::GridIndex;
469 using CoordScalar = typename GV::ctype;
470
471 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
472
473public:
475 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
477 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
479 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
480
482 using LocalView = typename Traits::template LocalView<ThisType, false>;
484 using SubControlVolume = typename Traits::SubControlVolume;
486 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
490 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
492 using DofMapper = typename Traits::ElementMapper;
494 using GridView = GV;
496 using MpfaHelper = typename Traits::template MpfaHelper<ThisType>;
497
500 static constexpr DiscretizationMethod discMethod{};
501
503 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
504
506 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
507
511 : ParentType(gridView)
512 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
513 { return is.boundary() || isBranching; } )
514 {
515 checkOverlapSizeCCMpfa(gridView);
516 update_();
517 }
518
520 CCMpfaFVGridGeometry(const GridView& gridView, const SecondaryIvIndicatorType& indicator)
521 : ParentType(gridView)
522 , secondaryIvIndicator_(indicator)
523 {
524 checkOverlapSizeCCMpfa(gridView);
525 update_();
526 }
527
530 const DofMapper& dofMapper() const
531 { return this->elementMapper(); }
532
534 std::size_t numScv() const
535 { return numScvs_; }
536
538 std::size_t numScvf() const
539 { return numScvf_; }
540
542 std::size_t numBoundaryScvf() const
543 { return numBoundaryScvf_; }
544
546 std::size_t numDofs() const
547 { return this->gridView().size(0); }
548
551 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary, bool> = 0>
552 bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
553 { return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
554
557 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary, bool> = 0>
558 constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
559 { return false; }
560
562 bool isGhostVertex(const Vertex& v) const
563 { return isGhostVertex_[this->vertexMapper().index(v)]; }
564
566 bool isGhostVertex(GridIndexType vIdxGlobal) const
567 { return isGhostVertex_[vIdxGlobal]; }
568
569
571 void update(const GridView& gridView)
572 {
573 ParentType::update(gridView);
574 update_();
575 }
576
578 void update(GridView&& gridView)
579 {
580 ParentType::update(std::move(gridView));
581 update_();
582 }
583
586 { return MpfaHelper(); }
587
589 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
590 { return scvfIndicesOfScv_[scvIdx]; }
591
593 const std::vector<ScvfOutsideGridIndexStorage>& neighborVolVarIndices(GridIndexType scvIdx) const
594 { return neighborVolVarIndices_[scvIdx]; }
595
598 const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
599 { return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
600
603 { return flipScvfIndices_; }
604
608 { return connectivityMap_; }
609
612 { return ivIndexSets_; }
613
614private:
615
616 void update_()
617 {
618
619 // stop the time required for the update
620 Dune::Timer timer;
621
622 // resize containers
623 numScvs_ = numDofs();
624 scvfIndicesOfScv_.resize(numScvs_);
625 neighborVolVarIndices_.resize(numScvs_);
626
627 // Some methods require to use a second type of interaction volume, e.g.
628 // around vertices on the boundary or branching points (surface grids)
629 const auto numVert = this->gridView().size(dim);
630 secondaryInteractionVolumeVertices_.assign(numVert, false);
631
632 // find vertices on processor boundaries HERE!!
633 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
634
635 // instantiate the dual grid index set (to be used for construction of interaction volumes)
636 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
637
638 // keep track of boundary scvfs and scvf vertex indices in order to set up flip scvf index set
639 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
640 std::vector<bool> scvfIsOnBoundary;
641 std::vector<GridIndexType> scvfVertexIndex;
642 scvfIsOnBoundary.reserve(maxNumScvfs);
643 scvfVertexIndex.reserve(maxNumScvfs);
644
645 // Build the SCVs and SCV faces
646 numScvf_ = 0;
647 numBoundaryScvf_ = 0;
648 for (const auto& element : elements(this->gridView()))
649 {
650 const auto eIdx = this->elementMapper().index(element);
651
652 // the element geometry
653 auto elementGeometry = element.geometry();
654
655 // the element-wise index sets for finite volume geometry
656 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
657 std::vector<GridIndexType> scvfsIndexSet;
658 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
659 scvfsIndexSet.reserve(numLocalFaces);
660 neighborVolVarIndexSet.reserve(numLocalFaces);
661
662 // for network grids there might be multiple intersections with the same geometryInInside
663 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
664 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
665 if (dim < dimWorld)
666 {
667 outsideIndices.resize(element.subEntities(1));
668 for (const auto& intersection : intersections(this->gridView(), element))
669 {
670 if (intersection.neighbor())
671 {
672 auto nIdx = this->elementMapper().index( intersection.outside() );
673 outsideIndices[intersection.indexInInside()].push_back(nIdx);
674 }
675 }
676 }
677
678 // construct the sub control volume faces
679 for (const auto& is : intersections(this->gridView(), element))
680 {
681 const auto indexInInside = is.indexInInside();
682 const bool boundary = is.boundary();
683 const bool neighbor = is.neighbor();
684
685 // for surface grids, skip the rest if handled already
686 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
687 continue;
688
689 // if outside level > inside level, use the outside element in the following
690 const bool useNeighbor = neighbor && is.outside().level() > element.level();
691 const auto& e = useNeighbor ? is.outside() : element;
692 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
693 const auto eg = e.geometry();
694 const auto refElement = referenceElement(eg);
695
696 // evaluate if vertices on this intersection use primary/secondary IVs
697 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
698 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
699
700 // make the scv faces belonging to each corner of the intersection
701 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
702 {
703 // get the global vertex index the scv face is connected to
704 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
705 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
706
707 // do not build scvfs connected to a processor boundary
708 if (isGhostVertex_[vIdxGlobal])
709 continue;
710
711 // if this vertex is tagged to use the secondary IVs, store info
712 if (usesSecondaryIV)
713 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
714
715 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
716 const auto& outsideScvIndices = [&] ()
717 {
718 if (!boundary)
719 return dim == dimWorld ?
720 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
721 outsideIndices[indexInInside];
722 else
723 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
724 } ();
725
726 // insert the scvf data into the dual grid index set
727 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
728
729 // store information on the scv face
730 scvfsIndexSet.push_back(numScvf_++);
731 scvfIsOnBoundary.push_back(boundary);
732 scvfVertexIndex.push_back(vIdxGlobal);
733 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
734 }
735
736 // for network grids, clear outside indices to not make a second scvf on that facet
737 if (dim < dimWorld)
738 outsideIndices[indexInInside].clear();
739 }
740
741 // store the sets of indices in the data container
742 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
743 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
744 }
745
746 // Make the flip scvf index set
747 flipScvfIndices_.resize(numScvf_);
748 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
749 {
750 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
751 for (unsigned int i = 0; i < scvfIndices.size(); ++i)
752 {
753 // boundary scvf have no flip scvfs
754 if (scvfIsOnBoundary[ scvfIndices[i] ])
755 continue;
756
757 const auto scvfIdx = scvfIndices[i];
758 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
759 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
760
761 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
762 for (unsigned int j = 0; j < numOutsideScvs; ++j)
763 {
764 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
765 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
766 for (unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
767 {
768 const auto outsideScvfIndex = outsideScvfIndices[k];
769 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
770 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
771 if (outsideScvfVertexIndex == vIdxGlobal &&
772 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
773 {
774 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
775 // there is always only one flip face in an outside element
776 break;
777 }
778 }
779 }
780 }
781 }
782
783 // building the geometries has finished
784 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
785
786 // Initialize the grid interaction volume index sets
787 timer.reset();
788 ivIndexSets_.update(*this, std::move(dualIdSet));
789 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
790
791 // build the connectivity map for an efficient assembly
792 timer.reset();
793 connectivityMap_.update(*this);
794 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
795 }
796
797 // connectivity map for efficient assembly
798 ConnectivityMap connectivityMap_;
799
800 // containers storing the global data
801 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
802 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
803 std::vector<bool> secondaryInteractionVolumeVertices_;
804 std::vector<bool> isGhostVertex_;
805 GridIndexType numScvs_;
806 GridIndexType numScvf_;
807 GridIndexType numBoundaryScvf_;
808
809 // needed for embedded surface and network grids (dim < dimWorld)
810 FlipScvfIndexSet flipScvfIndices_;
811
812 // The grid interaction volume index set
813 GridIVIndexSets ivIndexSets_;
814
815 // Indicator function on where to use the secondary IVs
816 SecondaryIvIndicatorType secondaryIvIndicator_;
817};
818
819} // end namespace Dumux
820
821#endif
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Check the overlap size for different discretization methods.
Helper classes to compute the integration elements.
Base class for grid geometries.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
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:232
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all grid geometries.
Definition: basegridgeometry.hh:61
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:69
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
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:222
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:213
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:205
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:178
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:185
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:200
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:209
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:217
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:226
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:192
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:196
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:458
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:479
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:486
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:490
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:566
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:598
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:578
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:611
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:488
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:477
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:593
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:520
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:602
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:510
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:534
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:589
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:475
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:484
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:542
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:492
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:585
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:530
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:607
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:562
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:496
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:482
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:552
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:558
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:538
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:571
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:546
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:40
Definition: method.hh:45