3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Loading...
Searching...
No Matches
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
35
36namespace Dumux {
37
49template<class GridView, class Traits, bool enableCache>
51
53template<class GridView>
54void checkOverlapSizeCCMpfa(const GridView& gridView)
55{
56 // Check if the overlap size is what we expect
58 DUNE_THROW(Dune::InvalidStateException, "The ccmpfa discretization method needs at least an overlap of 1 for parallel computations. "
59 << " Set the parameter \"Grid.Overlap\" in the input file.");
60}
61
68template<class GV, class Traits>
69class CCMpfaFVGridGeometry<GV, Traits, true>
70: public BaseGridGeometry<GV, Traits>
71{
73 using ParentType = BaseGridGeometry<GV, Traits>;
74
75 static constexpr int dim = GV::dimension;
76 static constexpr int dimWorld = GV::dimensionworld;
77
78 using Element = typename GV::template Codim<0>::Entity;
79 using Vertex = typename GV::template Codim<dim>::Entity;
80 using Intersection = typename GV::Intersection;
81 using GridIndexType = typename IndexTraits<GV>::GridIndex;
82 using CoordScalar = typename GV::ctype;
83
84 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
85
86public:
88 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
90 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
92 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
93
95 using LocalView = typename Traits::template LocalView<ThisType, true>;
97 using SubControlVolume = typename Traits::SubControlVolume;
99 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
132 : ParentType(gridView)
133 , secondaryIvIndicator_(indicator)
134 {
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 {
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 = referenceElement(eg);
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 = getParam<CoordScalar>("MPFA.Q");
282
283 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
284 const auto& outsideScvIndices = [&] ()
285 {
286 if (!boundary)
287 return dim == dimWorld ?
288 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
289 outsideIndices[indexInInside];
290 else
291 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs) + numBoundaryScvf_++});
292 } ();
293
294 scvfIndexSet.push_back(scvfIdx);
295 scvfs_.emplace_back(MpfaHelper(),
296 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
297 is,
298 vIdxGlobal,
299 vIdxLocal,
300 scvfIdx,
301 eIdx,
302 outsideScvIndices,
303 q,
304 boundary);
305
306 // insert the scvf data into the dual grid index set
307 dualIdSet[vIdxGlobal].insert(scvfs_.back());
308
309 // increment scvf counter
310 scvfIdx++;
311 }
312
313 // for network grids, clear outside indices to not make a second scvf on that facet
314 if (dim < dimWorld)
315 outsideIndices[indexInInside].clear();
316 }
317
318 // create sub control volume for this element
319 scvs_[eIdx] = SubControlVolume(std::move(elementGeometry), eIdx);
320
321 // Save the scvf indices belonging to this scv to build up fv element geometries fast
322 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
323 }
324
325 // Make the flip index set for network and surface grids
326 flipScvfIndices_.resize(scvfs_.size());
327 for (const auto& scvf : scvfs_)
328 {
329 if (scvf.boundary())
330 continue;
331
332 const auto numOutsideScvs = scvf.numOutsideScvs();
333 const auto vIdxGlobal = scvf.vertexIndex();
334 const auto insideScvIdx = scvf.insideScvIdx();
335
336 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
337 for (std::size_t i = 0; i < numOutsideScvs; ++i)
338 {
339 const auto outsideScvIdx = scvf.outsideScvIdx(i);
340 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
341 {
342 const auto& outsideScvf = this->scvf(outsideScvfIndex);
343 if (outsideScvf.vertexIndex() == vIdxGlobal &&
344 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
345 {
346 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
347 // there is always only one flip face in an outside element
348 break;
349 }
350 }
351 }
352 }
353
354 // building the geometries has finished
355 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
356
357 // Initialize the grid interaction volume index sets
358 timer.reset();
359 ivIndexSets_.update(*this, std::move(dualIdSet));
360 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
361
362 // build the connectivity map for an efficient assembly
363 timer.reset();
364 connectivityMap_.update(*this);
365 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
366 }
367
370 { return MpfaHelper(); }
371
373 const SubControlVolume& scv(GridIndexType scvIdx) const
374 { return scvs_[scvIdx]; }
375
377 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
378 { return scvfs_[scvfIdx]; }
379
383 { return connectivityMap_; }
384
387 { return ivIndexSets_; }
388
390 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
391 { return scvfIndicesOfScv_[scvIdx]; }
392
395 { return flipScvfIndices_; }
396
399 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
400 { return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]]; }
401
403 bool hasBoundaryScvf(GridIndexType eIdx) const
404 { return hasBoundaryScvf_[eIdx]; }
405
406private:
407 // connectivity map for efficient assembly
408 ConnectivityMap connectivityMap_;
409
410 // the finite volume grid geometries
411 std::vector<SubControlVolume> scvs_;
412 std::vector<SubControlVolumeFace> scvfs_;
413
414 // containers storing the global data
415 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
416 std::vector<bool> secondaryInteractionVolumeVertices_;
417 GridIndexType numBoundaryScvf_;
418 std::vector<bool> hasBoundaryScvf_;
419
420 // needed for embedded surface and network grids (dim < dimWorld)
421 FlipScvfIndexSet flipScvfIndices_;
422
423 // The grid interaction volume index set
424 GridIVIndexSets ivIndexSets_;
425
426 // Indicator function on where to use the secondary IVs
427 SecondaryIvIndicatorType secondaryIvIndicator_;
428};
429
437template<class GV, class Traits>
438class CCMpfaFVGridGeometry<GV, Traits, false>
439: public BaseGridGeometry<GV, Traits>
440{
442 using ParentType = BaseGridGeometry<GV, Traits>;
443
444 static constexpr int dim = GV::dimension;
445 static constexpr int dimWorld = GV::dimensionworld;
446
447 using Element = typename GV::template Codim<0>::Entity;
448 using Vertex = typename GV::template Codim<dim>::Entity;
449 using Intersection = typename GV::Intersection;
450 using GridIndexType = typename IndexTraits<GV>::GridIndex;
451 using CoordScalar = typename GV::ctype;
452
453 using ScvfOutsideGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
454
455public:
457 using FlipScvfIndexSet = std::vector<ScvfOutsideGridIndexStorage>;
459 using GridIVIndexSets = typename Traits::template GridIvIndexSets<ThisType>;
461 using SecondaryIvIndicatorType = std::function<bool(const Element&, const Intersection&, bool)>;
462
464 using LocalView = typename Traits::template LocalView<ThisType, false>;
466 using SubControlVolume = typename Traits::SubControlVolume;
468 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
472 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
474 using DofMapper = typename Traits::ElementMapper;
476 using GridView = GV;
478 using MpfaHelper = typename Traits::template MpfaHelper<ThisType>;
479
482
484 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
485
487 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
488
492 : ParentType(gridView)
493 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
494 { return is.boundary() || isBranching; } )
495 {
496 checkOverlapSizeCCMpfa(gridView);
497 }
498
501 : ParentType(gridView)
502 , secondaryIvIndicator_(indicator)
503 {
505 }
506
509 const DofMapper& dofMapper() const
510 { return this->elementMapper(); }
511
513 std::size_t numScv() const
514 { return numScvs_; }
515
517 std::size_t numScvf() const
518 { return numScvf_; }
519
521 std::size_t numBoundaryScvf() const
522 { return numBoundaryScvf_; }
523
525 std::size_t numDofs() const
526 { return this->gridView().size(0); }
527
530 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary, bool> = 0>
531 bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
532 { return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
533
536 template<bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary, bool> = 0>
537 constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
538 { return false; }
539
541 bool isGhostVertex(const Vertex& v) const
542 { return isGhostVertex_[this->vertexMapper().index(v)]; }
543
545 bool isGhostVertex(GridIndexType vIdxGlobal) const
546 { return isGhostVertex_[vIdxGlobal]; }
547
549 void update()
550 {
552
553 // stop the time required for the update
554 Dune::Timer timer;
555
556 // resize containers
557 numScvs_ = numDofs();
558 scvfIndicesOfScv_.resize(numScvs_);
559 neighborVolVarIndices_.resize(numScvs_);
560
561 // Some methods require to use a second type of interaction volume, e.g.
562 // around vertices on the boundary or branching points (surface grids)
563 const auto numVert = this->gridView().size(dim);
564 secondaryInteractionVolumeVertices_.resize(numVert, false);
565
566 // find vertices on processor boundaries HERE!!
567 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
568
569 // instantiate the dual grid index set (to be used for construction of interaction volumes)
570 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
571
572 // keep track of boundary scvfs and scvf vertex indices in order to set up flip scvf index set
573 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
574 std::vector<bool> scvfIsOnBoundary;
575 std::vector<GridIndexType> scvfVertexIndex;
576 scvfIsOnBoundary.reserve(maxNumScvfs);
577 scvfVertexIndex.reserve(maxNumScvfs);
578
579 // Build the SCVs and SCV faces
580 numScvf_ = 0;
581 numBoundaryScvf_ = 0;
582 for (const auto& element : elements(this->gridView()))
583 {
584 const auto eIdx = this->elementMapper().index(element);
585
586 // the element geometry
587 auto elementGeometry = element.geometry();
588
589 // the element-wise index sets for finite volume geometry
590 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
591 std::vector<GridIndexType> scvfsIndexSet;
592 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
593 scvfsIndexSet.reserve(numLocalFaces);
594 neighborVolVarIndexSet.reserve(numLocalFaces);
595
596 // for network grids there might be multiple intersections with the same geometryInInside
597 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
598 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
599 if (dim < dimWorld)
600 {
601 outsideIndices.resize(element.subEntities(1));
602 for (const auto& intersection : intersections(this->gridView(), element))
603 {
604 if (intersection.neighbor())
605 {
606 auto nIdx = this->elementMapper().index( intersection.outside() );
607 outsideIndices[intersection.indexInInside()].push_back(nIdx);
608 }
609 }
610 }
611
612 // construct the sub control volume faces
613 for (const auto& is : intersections(this->gridView(), element))
614 {
615 const auto indexInInside = is.indexInInside();
616 const bool boundary = is.boundary();
617 const bool neighbor = is.neighbor();
618
619 // for surface grids, skip the rest if handled already
620 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
621 continue;
622
623 // if outside level > inside level, use the outside element in the following
624 const bool useNeighbor = neighbor && is.outside().level() > element.level();
625 const auto& e = useNeighbor ? is.outside() : element;
626 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
627 const auto eg = e.geometry();
628 const auto refElement = referenceElement(eg);
629
630 // evaluate if vertices on this intersection use primary/secondary IVs
631 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 : false;
632 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
633
634 // make the scv faces belonging to each corner of the intersection
635 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
636 {
637 // get the global vertex index the scv face is connected to
638 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
639 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
640
641 // do not build scvfs connected to a processor boundary
642 if (isGhostVertex_[vIdxGlobal])
643 continue;
644
645 // if this vertex is tagged to use the secondary IVs, store info
646 if (usesSecondaryIV)
647 secondaryInteractionVolumeVertices_[vIdxGlobal] = true;
648
649 // make the scv face (for non-boundary scvfs on network grids, use precalculated outside indices)
650 const auto& outsideScvIndices = [&] ()
651 {
652 if (!boundary)
653 return dim == dimWorld ?
654 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
655 outsideIndices[indexInInside];
656 else
657 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
658 } ();
659
660 // insert the scvf data into the dual grid index set
661 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
662
663 // store information on the scv face
664 scvfsIndexSet.push_back(numScvf_++);
665 scvfIsOnBoundary.push_back(boundary);
666 scvfVertexIndex.push_back(vIdxGlobal);
667 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
668 }
669
670 // for network grids, clear outside indices to not make a second scvf on that facet
671 if (dim < dimWorld)
672 outsideIndices[indexInInside].clear();
673 }
674
675 // store the sets of indices in the data container
676 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
677 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
678 }
679
680 // Make the flip scvf index set
681 flipScvfIndices_.resize(numScvf_);
682 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
683 {
684 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
685 for (unsigned int i = 0; i < scvfIndices.size(); ++i)
686 {
687 // boundary scvf have no flip scvfs
688 if (scvfIsOnBoundary[ scvfIndices[i] ])
689 continue;
690
691 const auto scvfIdx = scvfIndices[i];
692 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
693 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
694
695 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
696 for (unsigned int j = 0; j < numOutsideScvs; ++j)
697 {
698 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
699 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
700 for (unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
701 {
702 const auto outsideScvfIndex = outsideScvfIndices[k];
703 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
704 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
705 if (outsideScvfVertexIndex == vIdxGlobal &&
706 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
707 {
708 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
709 // there is always only one flip face in an outside element
710 break;
711 }
712 }
713 }
714 }
715 }
716
717 // building the geometries has finished
718 std::cout << "Initializing of the grid finite volume geometry took " << timer.elapsed() << " seconds." << std::endl;
719
720 // Initialize the grid interaction volume index sets
721 timer.reset();
722 ivIndexSets_.update(*this, std::move(dualIdSet));
723 std::cout << "Initializing of the grid interaction volume index sets took " << timer.elapsed() << " seconds." << std::endl;
724
725 // build the connectivity map for an efficient assembly
726 timer.reset();
727 connectivityMap_.update(*this);
728 std::cout << "Initializing of the connectivity map took " << timer.elapsed() << " seconds." << std::endl;
729 }
730
733 { return MpfaHelper(); }
734
736 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
737 { return scvfIndicesOfScv_[scvIdx]; }
738
740 const std::vector<ScvfOutsideGridIndexStorage>& neighborVolVarIndices(GridIndexType scvIdx) const
741 { return neighborVolVarIndices_[scvIdx]; }
742
745 const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
746 { return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
747
750 { return flipScvfIndices_; }
751
755 { return connectivityMap_; }
756
759 { return ivIndexSets_; }
760
761private:
762 // connectivity map for efficient assembly
763 ConnectivityMap connectivityMap_;
764
765 // containers storing the global data
766 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
767 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
768 std::vector<bool> secondaryInteractionVolumeVertices_;
769 std::vector<bool> isGhostVertex_;
770 GridIndexType numScvs_;
771 GridIndexType numScvf_;
772 GridIndexType numBoundaryScvf_;
773
774 // needed for embedded surface and network grids (dim < dimWorld)
775 FlipScvfIndexSet flipScvfIndices_;
776
777 // The grid interaction volume index set
778 GridIVIndexSets ivIndexSets_;
779
780 // Indicator function on where to use the secondary IVs
781 SecondaryIvIndicatorType secondaryIvIndicator_;
782};
783
784} // end namespace Dumux
785
786#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.
Helper classes to compute the integration elements.
Base class for grid geometries.
BaseGridGeometry(const GridView &gridView)
Constructor computes the bouding box of the entire domain, for e.g. setting boundary conditions.
Definition basegridgeometry.hh:78
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ ccmpfa
Definition method.hh:38
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:348
Definition adapt.hh:29
void checkOverlapSizeCCMpfa(const GridView &gridView)
check the overlap size for parallel computations
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:54
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:177
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:39
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices for constant grids.
Definition basegridgeometry.hh:127
Element element(GridIndexType eIdx) const
Get an element from a global element index.
Definition basegridgeometry.hh:171
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices for constant grids.
Definition basegridgeometry.hh:121
const GridView & gridView() const
Return the gridView this grid geometry object lives on.
Definition basegridgeometry.hh:115
void update()
Update all fvElementGeometries (do this again after grid adaption).
Definition basegridgeometry.hh:91
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:50
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:105
static constexpr int maxElementStencilSize
The maximum admissible stencil size (used for static memory allocation during assembly).
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:115
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:99
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:399
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:390
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:101
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:162
static constexpr bool hasSingleInteractionVolumeType
State if only a single type is used for interaction volumes.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:118
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:88
const ConnectivityMap & connectivityMap() const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:382
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:97
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:90
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:95
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:377
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:386
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:394
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:92
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:103
GV GridView
export the grid view type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:107
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:403
CCMpfaFVGridGeometry(const GridView &gridView)
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:122
static constexpr DiscretizationMethod discMethod
export the discretization method this geometry belongs to
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:112
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:369
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:373
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:168
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:461
static constexpr bool hasSingleInteractionVolumeType
State if only a single type is used for interaction volumes.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:487
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:468
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:472
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:545
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:745
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:758
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:470
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:459
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:740
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:500
static constexpr DiscretizationMethod discMethod
export the discretization method this geometry belongs to
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:481
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:749
static constexpr int maxElementStencilSize
The maximum admissible stencil size (used for static memory allocation during assembly).
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:484
CCMpfaFVGridGeometry(const GridView &gridView)
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:491
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:513
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:736
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:457
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:466
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:521
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:474
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:732
GV GridView
export the grid view type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:476
void update()
Updates all finite volume geometries of the grid. Has to be called again after grid adaption.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:549
const DofMapper & dofMapper() const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:509
const ConnectivityMap & connectivityMap() const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:754
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:541
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:478
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:464
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:531
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:537
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:517
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition discretization/cellcentered/mpfa/fvgridgeometry.hh:525
static bool isValid(const GridView &gridView) noexcept
Definition checkoverlapsize.hh:42