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