3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/cellcentered/tpfa/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_CCTPFA_FV_GRID_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCTPFA_FV_GRID_GEOMETRY_HH
28
29#include <algorithm>
30
33
42
43namespace Dumux {
44
51template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
53: public MapperTraits
54{
57
58 template<class GridGeometry>
60
61 template<class GridGeometry, bool enableCache>
63
68 static constexpr int maxNumScvfNeighbors = int(GridView::dimension)<int(GridView::dimensionworld) ? 8 : 1<<(GridView::dimension-1);
69};
70
77template<class GridView,
78 bool enableGridGeometryCache = false,
81
88template<class GV, class Traits>
89class CCTpfaFVGridGeometry<GV, true, Traits>
90: public BaseGridGeometry<GV, Traits>
91{
94 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
95 using GridIndexType = typename IndexTraits<GV>::GridIndex;
96 using Element = typename GV::template Codim<0>::Entity;
97
98 static const int dim = GV::dimension;
99 static const int dimWorld = GV::dimensionworld;
100
101public:
103 using LocalView = typename Traits::template LocalView<ThisType, true>;
105 using SubControlVolume = typename Traits::SubControlVolume;
107 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
111 using DofMapper = typename Traits::ElementMapper;
112
115
117 static constexpr int maxElementStencilSize = LocalView::maxNumElementScvfs*Traits::maxNumScvfNeighbors + 1;
118
120 using GridView = GV;
121
124 : ParentType(gridView)
125 {
126 // Check if the overlap size is what we expect
128 DUNE_THROW(Dune::InvalidStateException, "The cctpfa discretization method needs at least an overlap of 1 for parallel computations. "
129 << " Set the parameter \"Grid.Overlap\" in the input file.");
130 }
131
134 const DofMapper& dofMapper() const
135 { return this->elementMapper(); }
136
138 std::size_t numScv() const
139 {
140 return scvs_.size();
141 }
142
144 std::size_t numScvf() const
145 {
146 return scvfs_.size();
147 }
148
150 std::size_t numBoundaryScvf() const
151 {
152 return numBoundaryScvf_;
153 }
154
156 std::size_t numDofs() const
157 { return this->gridView().size(0); }
158
160 void update()
161 {
162 ParentType::update();
163
164 // clear containers (necessary after grid refinement)
165 scvs_.clear();
166 scvfs_.clear();
167 scvfIndicesOfScv_.clear();
168 flipScvfIndices_.clear();
169
170 // determine size of containers
171 std::size_t numScvs = numDofs();
172 std::size_t numScvf = 0;
173 for (const auto& element : elements(this->gridView()))
174 numScvf += element.subEntities(1);
175
176 // reserve memory
177 scvs_.resize(numScvs);
178 scvfs_.reserve(numScvf);
179 scvfIndicesOfScv_.resize(numScvs);
180 hasBoundaryScvf_.resize(numScvs, false);
181
182 // Build the scvs and scv faces
183 GridIndexType scvfIdx = 0;
184 numBoundaryScvf_ = 0;
185 for (const auto& element : elements(this->gridView()))
186 {
187 const auto eIdx = this->elementMapper().index(element);
188 scvs_[eIdx] = SubControlVolume(element.geometry(), eIdx);
189
190 // the element-wise index sets for finite volume geometry
191 std::vector<GridIndexType> scvfsIndexSet;
192 scvfsIndexSet.reserve(element.subEntities(1));
193
194 // for network grids there might be multiple intersection with the same geometryInInside
195 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
196 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
197 std::vector<ScvfGridIndexStorage> outsideIndices;
198 if (dim < dimWorld)
199 {
201 outsideIndices.resize(element.subEntities(1));
202 std::for_each(outsideIndices.begin(), outsideIndices.end(), [eIdx] (auto& nIndices) { nIndices.push_back(eIdx); });
203
204 // second, insert neighbors
205 for (const auto& intersection : intersections(this->gridView(), element))
206 {
207 if (intersection.neighbor())
208 {
209 const auto nIdx = this->elementMapper().index( intersection.outside() );
210 outsideIndices[intersection.indexInInside()].push_back(nIdx);
211 }
212 }
213 }
214
215 for (const auto& intersection : intersections(this->gridView(), element))
216 {
217 // inner sub control volume faces (includes periodic boundaries)
218 if (intersection.neighbor())
219 {
220 // update the grid geometry if we have periodic boundaries
221 if (intersection.boundary())
222 this->setPeriodic();
223
224 if (dim == dimWorld)
225 {
226 const auto nIdx = this->elementMapper().index(intersection.outside());
227 scvfs_.emplace_back(intersection,
228 intersection.geometry(),
229 scvfIdx,
230 ScvfGridIndexStorage({eIdx, nIdx}),
231 false);
232 scvfsIndexSet.push_back(scvfIdx++);
233 }
234 // this is for network grids
235 // (will be optimized away of dim == dimWorld)
236 else
237 {
238 auto indexInInside = intersection.indexInInside();
239 // check if we already handled this facet
240 if (outsideIndices[indexInInside].empty())
241 continue;
242 else
243 {
244 scvfs_.emplace_back(intersection,
245 intersection.geometry(),
246 scvfIdx,
247 outsideIndices[indexInInside],
248 false);
249 scvfsIndexSet.push_back(scvfIdx++);
250 outsideIndices[indexInInside].clear();
251 }
252 }
253 }
254 // boundary sub control volume faces
255 else if (intersection.boundary())
256 {
257 scvfs_.emplace_back(intersection,
258 intersection.geometry(),
259 scvfIdx,
260 ScvfGridIndexStorage({eIdx, static_cast<GridIndexType>(this->gridView().size(0) + numBoundaryScvf_++)}),
261 true);
262 scvfsIndexSet.push_back(scvfIdx++);
263
264 hasBoundaryScvf_[eIdx] = true;
265 }
266 }
267
268 // Save the scvf indices belonging to this scv to build up fv element geometries fast
269 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
270 }
271
272 // Make the flip index set for network, surface, and periodic grids
273 if (dim < dimWorld || this->isPeriodic())
274 {
275 flipScvfIndices_.resize(scvfs_.size());
276 for (auto&& scvf : scvfs_)
277 {
278 if (scvf.boundary())
279 continue;
280
281 flipScvfIndices_[scvf.index()].resize(scvf.numOutsideScvs());
282 const auto insideScvIdx = scvf.insideScvIdx();
283 // check which outside scvf has the insideScvIdx index in its outsideScvIndices
284 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
285 flipScvfIndices_[scvf.index()][i] = findFlippedScvfIndex_(insideScvIdx, scvf.outsideScvIdx(i));
286 }
287 }
288
289 // build the connectivity map for an effecient assembly
290 connectivityMap_.update(*this);
291 }
292
294 const SubControlVolume& scv(GridIndexType scvIdx) const
295 {
296 return scvs_[scvIdx];
297 }
298
300 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
301 {
302 return scvfs_[scvfIdx];
303 }
304
307 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
308 {
309 return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]];
310 }
311
313 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
314 {
315 return scvfIndicesOfScv_[scvIdx];
316 }
317
322 const ConnectivityMap &connectivityMap() const
323 { return connectivityMap_; }
324
326 bool hasBoundaryScvf(GridIndexType eIdx) const
327 { return hasBoundaryScvf_[eIdx]; }
328
329private:
330 // find the scvf that has insideScvIdx in its outsideScvIdx list and outsideScvIdx as its insideScvIdx
331 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType outsideScvIdx)
332 {
333 // go over all potential scvfs of the outside scv
334 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
335 {
336 const auto& outsideScvf = this->scvf(outsideScvfIndex);
337 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
338 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
339 return outsideScvf.index();
340 }
341
342 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
343 }
344
346 ConnectivityMap connectivityMap_;
347
349 std::vector<SubControlVolume> scvs_;
350 std::vector<SubControlVolumeFace> scvfs_;
351 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
352 std::size_t numBoundaryScvf_;
353 std::vector<bool> hasBoundaryScvf_;
354
356 std::vector<std::vector<GridIndexType>> flipScvfIndices_;
357};
358
366template<class GV, class Traits>
367class CCTpfaFVGridGeometry<GV, false, Traits>
368: public BaseGridGeometry<GV, Traits>
369{
372 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
373
374 using GridIndexType = typename IndexTraits<GV>::GridIndex;
375 using Element = typename GV::template Codim<0>::Entity;
376
377 static const int dim = GV::dimension;
378 static const int dimWorld = GV::dimensionworld;
379
380 using ScvfGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::GridIndexStorage;
381 using NeighborVolVarIndices = typename std::conditional_t< (dim<dimWorld),
382 ScvfGridIndexStorage,
383 Dune::ReservedVector<GridIndexType, 1> >;
384
385public:
387 using LocalView = typename Traits::template LocalView<ThisType, false>;
389 using SubControlVolume = typename Traits::SubControlVolume;
391 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
395 using DofMapper = typename Traits::ElementMapper;
396
399
401 static constexpr int maxElementStencilSize = LocalView::maxNumElementScvfs*Traits::maxNumScvfNeighbors + 1;
402
404 using GridView = GV;
405
408 : ParentType(gridView)
409 {
410 // Check if the overlap size is what we expect
412 DUNE_THROW(Dune::InvalidStateException, "The cctpfa discretization method needs at least an overlap of 1 for parallel computations. "
413 << " Set the parameter \"Grid.Overlap\" in the input file.");
414 }
415
418 const DofMapper& dofMapper() const
419 { return this->elementMapper(); }
420
422 std::size_t numScv() const
423 {
424 return numScvs_;
425 }
426
428 std::size_t numScvf() const
429 {
430 return numScvf_;
431 }
432
434 std::size_t numBoundaryScvf() const
435 {
436 return numBoundaryScvf_;
437 }
438
440 std::size_t numDofs() const
441 { return this->gridView().size(0); }
442
444 void update()
445 {
446 ParentType::update();
447
448 // clear local data
449 scvfIndicesOfScv_.clear();
450 neighborVolVarIndices_.clear();
451
452 // reserve memory or resize the containers
453 numScvs_ = numDofs();
454 numScvf_ = 0;
455 numBoundaryScvf_ = 0;
456 scvfIndicesOfScv_.resize(numScvs_);
457 neighborVolVarIndices_.resize(numScvs_);
458
459 // Build the SCV and SCV face
460 for (const auto& element : elements(this->gridView()))
461 {
462 const auto eIdx = this->elementMapper().index(element);
463
464 // the element-wise index sets for finite volume geometry
465 auto numLocalFaces = element.subEntities(1);
466 std::vector<GridIndexType> scvfsIndexSet;
467 std::vector<NeighborVolVarIndices> neighborVolVarIndexSet;
468 scvfsIndexSet.reserve(numLocalFaces);
469 neighborVolVarIndexSet.reserve(numLocalFaces);
470
471 // for network grids there might be multiple intersection with the same geometryInInside
472 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
473 std::vector<NeighborVolVarIndices> outsideIndices;
474 if (dim < dimWorld)
475 {
476 outsideIndices.resize(numLocalFaces);
477 for (const auto& intersection : intersections(this->gridView(), element))
478 {
479 if (intersection.neighbor())
480 {
481 const auto nIdx = this->elementMapper().index(intersection.outside());
482 outsideIndices[intersection.indexInInside()].push_back(nIdx);
483 }
484 }
485 }
486
487 for (const auto& intersection : intersections(this->gridView(), element))
488 {
489 // inner sub control volume faces (includes periodic boundaries)
490 if (intersection.neighbor())
491 {
492 // update the grid geometry if we have periodic boundaries
493 if (intersection.boundary())
494 this->setPeriodic();
495
496 if (dim == dimWorld)
497 {
498 scvfsIndexSet.push_back(numScvf_++);
499 const auto nIdx = this->elementMapper().index(intersection.outside());
500 neighborVolVarIndexSet.emplace_back(NeighborVolVarIndices({nIdx}));
501 }
502 // this is for network grids
503 // (will be optimized away of dim == dimWorld)
504 else
505 {
506 auto indexInInside = intersection.indexInInside();
507 // check if we already handled this facet
508 if (outsideIndices[indexInInside].empty())
509 continue;
510 else
511 {
512 scvfsIndexSet.push_back(numScvf_++);
513 neighborVolVarIndexSet.emplace_back(std::move(outsideIndices[indexInInside]));
514 outsideIndices[indexInInside].clear();
515 }
516 }
517 }
518 // boundary sub control volume faces
519 else if (intersection.boundary())
520 {
521 scvfsIndexSet.push_back(numScvf_++);
522 neighborVolVarIndexSet.emplace_back(NeighborVolVarIndices({static_cast<GridIndexType>(numScvs_ + numBoundaryScvf_++)}));
523 }
524 }
525
526 // store the sets of indices in the data container
527 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
528 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
529 }
530
531 // build the connectivity map for an effecient assembly
532 connectivityMap_.update(*this);
533 }
534
535 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
536 { return scvfIndicesOfScv_[scvIdx]; }
537
539 const std::vector<NeighborVolVarIndices>& neighborVolVarIndices(GridIndexType scvIdx) const
540 { return neighborVolVarIndices_[scvIdx]; }
541
546 const ConnectivityMap &connectivityMap() const
547 { return connectivityMap_; }
548
549private:
550
552 std::size_t numScvs_;
553 std::size_t numScvf_;
554 std::size_t numBoundaryScvf_;
555
557 ConnectivityMap connectivityMap_;
558
560 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
561 std::vector<std::vector<NeighborVolVarIndices>> neighborVolVarIndices_;
562};
563
564} // end namespace Dumux
565
566#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
Check the overlap size for different discretization methods.
The available discretization methods in Dumux.
Base class for grid geometries.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
A simple version of the connectivity map for cellcentered schemes. This implementation works for sche...
Definition: cellcentered/connectivitymap.hh:53
Sub control volumes for cell-centered discretization schemes.
Definition: discretization/cellcentered/subcontrolvolume.hh:63
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:50
The default traits for the tpfa finite volume grid geometry Defines the scv and scvf types and the ma...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:54
static constexpr int maxNumScvfNeighbors
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:68
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:80
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:91
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:103
const ConnectivityMap & connectivityMap() const
Returns the connectivity map of which dofs have derivatives with respect to a given dof.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:322
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:160
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:138
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:109
CCTpfaFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:123
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:105
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:107
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:156
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Get the sub control volume face indices of an scv by global index.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:313
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:144
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:300
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:326
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:150
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:307
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:111
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:294
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:134
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:369
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:387
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:535
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:395
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:418
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:391
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:389
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:428
const ConnectivityMap & connectivityMap() const
Returns the connectivity map of which dofs have derivatives with respect to a given dof.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:546
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:444
CCTpfaFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:407
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:393
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:440
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:434
const std::vector< NeighborVolVarIndices > & neighborVolVarIndices(GridIndexType scvIdx) const
Return the neighbor volVar indices for all scvfs in the scv with index scvIdx.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:539
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:422
The sub control volume face.
Definition: discretization/cellcentered/tpfa/subcontrolvolumeface.hh:90
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:40
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Sub control volumes for cell-centered discretization schemes.
The sub control volume face.
Stores the face indices corresponding to the neighbors of an element that contribute to the derivativ...