3.5-git
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 <utility>
30#include <algorithm>
31
34
43
44namespace Dumux {
45
52template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
54: public MapperTraits
55{
58
59 template<class GridGeometry>
61
62 template<class GridGeometry, bool enableCache>
64
69 static constexpr int maxNumScvfNeighbors = int(GridView::dimension)<int(GridView::dimensionworld) ? 8 : 1<<(GridView::dimension-1);
70};
71
78template<class GridView,
79 bool enableGridGeometryCache = false,
82
89template<class GV, class Traits>
90class CCTpfaFVGridGeometry<GV, true, Traits>
91: public BaseGridGeometry<GV, Traits>
92{
95 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
96 using GridIndexType = typename IndexTraits<GV>::GridIndex;
97 using Element = typename GV::template Codim<0>::Entity;
98
99 static const int dim = GV::dimension;
100 static const int dimWorld = GV::dimensionworld;
101
102public:
104 using LocalView = typename Traits::template LocalView<ThisType, true>;
106 using SubControlVolume = typename Traits::SubControlVolume;
108 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
112 using DofMapper = typename Traits::ElementMapper;
113
116 static constexpr DiscretizationMethod discMethod{};
117
119 static constexpr int maxElementStencilSize = LocalView::maxNumElementScvfs*Traits::maxNumScvfNeighbors + 1;
120
122 using GridView = GV;
123
126 : ParentType(gridView)
127 {
128 // Check if the overlap size is what we expect
130 DUNE_THROW(Dune::InvalidStateException, "The cctpfa discretization method needs at least an overlap of 1 for parallel computations. "
131 << " Set the parameter \"Grid.Overlap\" in the input file.");
132
133 update_();
134 }
135
138 const DofMapper& dofMapper() const
139 { return this->elementMapper(); }
140
142 std::size_t numScv() const
143 {
144 return scvs_.size();
145 }
146
148 std::size_t numScvf() const
149 {
150 return scvfs_.size();
151 }
152
154 std::size_t numBoundaryScvf() const
155 {
156 return numBoundaryScvf_;
157 }
158
160 std::size_t numDofs() const
161 { return this->gridView().size(0); }
162
164 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
165 void update()
166 {
167 ParentType::update();
168 update_();
169 }
170
172 void update(const GridView& gridView)
173 {
174 ParentType::update(gridView);
175 update_();
176 }
177
179 void update(GridView&& gridView)
180 {
181 ParentType::update(std::move(gridView));
182 update_();
183 }
184
186 const SubControlVolume& scv(GridIndexType scvIdx) const
187 {
188 return scvs_[scvIdx];
189 }
190
192 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
193 {
194 return scvfs_[scvfIdx];
195 }
196
199 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx = 0) const
200 {
201 return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]];
202 }
203
205 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
206 {
207 return scvfIndicesOfScv_[scvIdx];
208 }
209
214 const ConnectivityMap &connectivityMap() const
215 { return connectivityMap_; }
216
218 bool hasBoundaryScvf(GridIndexType eIdx) const
219 { return hasBoundaryScvf_[eIdx]; }
220
221private:
222
223 void update_()
224 {
225 // clear containers (necessary after grid refinement)
226 scvs_.clear();
227 scvfs_.clear();
228 scvfIndicesOfScv_.clear();
229 flipScvfIndices_.clear();
230
231 // determine size of containers
232 std::size_t numScvs = numDofs();
233 std::size_t numScvf = 0;
234 for (const auto& element : elements(this->gridView()))
235 numScvf += element.subEntities(1);
236
237 // reserve memory
238 scvs_.resize(numScvs);
239 scvfs_.reserve(numScvf);
240 scvfIndicesOfScv_.resize(numScvs);
241 hasBoundaryScvf_.assign(numScvs, false);
242
243 // Build the scvs and scv faces
244 GridIndexType scvfIdx = 0;
245 numBoundaryScvf_ = 0;
246 for (const auto& element : elements(this->gridView()))
247 {
248 const auto eIdx = this->elementMapper().index(element);
249 scvs_[eIdx] = SubControlVolume(element.geometry(), eIdx);
250
251 // the element-wise index sets for finite volume geometry
252 std::vector<GridIndexType> scvfsIndexSet;
253 scvfsIndexSet.reserve(element.subEntities(1));
254
255 // for network grids there might be multiple intersection with the same geometryInInside
256 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
257 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
258 std::vector<ScvfGridIndexStorage> outsideIndices;
259 if (dim < dimWorld)
260 {
262 outsideIndices.resize(element.subEntities(1));
263 std::for_each(outsideIndices.begin(), outsideIndices.end(), [eIdx] (auto& nIndices) { nIndices.push_back(eIdx); });
264
265 // second, insert neighbors
266 for (const auto& intersection : intersections(this->gridView(), element))
267 {
268 if (intersection.neighbor())
269 {
270 const auto nIdx = this->elementMapper().index( intersection.outside() );
271 outsideIndices[intersection.indexInInside()].push_back(nIdx);
272 }
273 }
274 }
275
276 for (const auto& intersection : intersections(this->gridView(), element))
277 {
278 // inner sub control volume faces (includes periodic boundaries)
279 if (intersection.neighbor())
280 {
281 // update the grid geometry if we have periodic boundaries
282 if (intersection.boundary())
283 this->setPeriodic();
284
285 if (dim == dimWorld)
286 {
287 const auto nIdx = this->elementMapper().index(intersection.outside());
288 scvfs_.emplace_back(intersection,
289 intersection.geometry(),
290 scvfIdx,
291 ScvfGridIndexStorage({eIdx, nIdx}),
292 false);
293 scvfsIndexSet.push_back(scvfIdx++);
294 }
295 // this is for network grids
296 // (will be optimized away of dim == dimWorld)
297 else
298 {
299 auto indexInInside = intersection.indexInInside();
300 // check if we already handled this facet
301 if (outsideIndices[indexInInside].empty())
302 continue;
303 else
304 {
305 scvfs_.emplace_back(intersection,
306 intersection.geometry(),
307 scvfIdx,
308 outsideIndices[indexInInside],
309 false);
310 scvfsIndexSet.push_back(scvfIdx++);
311 outsideIndices[indexInInside].clear();
312 }
313 }
314 }
315 // boundary sub control volume faces
316 else if (intersection.boundary())
317 {
318 scvfs_.emplace_back(intersection,
319 intersection.geometry(),
320 scvfIdx,
321 ScvfGridIndexStorage({eIdx, static_cast<GridIndexType>(this->gridView().size(0) + numBoundaryScvf_++)}),
322 true);
323 scvfsIndexSet.push_back(scvfIdx++);
324
325 hasBoundaryScvf_[eIdx] = true;
326 }
327 }
328
329 // Save the scvf indices belonging to this scv to build up fv element geometries fast
330 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
331 }
332
333 // Make the flip index set for network, surface, and periodic grids
334 if (dim < dimWorld || this->isPeriodic())
335 {
336 flipScvfIndices_.resize(scvfs_.size());
337 for (auto&& scvf : scvfs_)
338 {
339 if (scvf.boundary())
340 continue;
341
342 flipScvfIndices_[scvf.index()].resize(scvf.numOutsideScvs());
343 const auto insideScvIdx = scvf.insideScvIdx();
344 // check which outside scvf has the insideScvIdx index in its outsideScvIndices
345 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
346 flipScvfIndices_[scvf.index()][i] = findFlippedScvfIndex_(insideScvIdx, scvf.outsideScvIdx(i));
347 }
348 }
349
350 // build the connectivity map for an efficient assembly
351 connectivityMap_.update(*this);
352 }
353
354 // find the scvf that has insideScvIdx in its outsideScvIdx list and outsideScvIdx as its insideScvIdx
355 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType outsideScvIdx)
356 {
357 // go over all potential scvfs of the outside scv
358 for (auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
359 {
360 const auto& outsideScvf = this->scvf(outsideScvfIndex);
361 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
362 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
363 return outsideScvf.index();
364 }
365
366 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
367 }
368
370 ConnectivityMap connectivityMap_;
371
373 std::vector<SubControlVolume> scvs_;
374 std::vector<SubControlVolumeFace> scvfs_;
375 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
376 std::size_t numBoundaryScvf_;
377 std::vector<bool> hasBoundaryScvf_;
378
380 std::vector<std::vector<GridIndexType>> flipScvfIndices_;
381};
382
390template<class GV, class Traits>
391class CCTpfaFVGridGeometry<GV, false, Traits>
392: public BaseGridGeometry<GV, Traits>
393{
396 using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>;
397
398 using GridIndexType = typename IndexTraits<GV>::GridIndex;
399 using Element = typename GV::template Codim<0>::Entity;
400
401 static const int dim = GV::dimension;
402 static const int dimWorld = GV::dimensionworld;
403
404 using ScvfGridIndexStorage = typename Traits::SubControlVolumeFace::Traits::GridIndexStorage;
405 using NeighborVolVarIndices = typename std::conditional_t< (dim<dimWorld),
406 ScvfGridIndexStorage,
407 Dune::ReservedVector<GridIndexType, 1> >;
408
409public:
411 using LocalView = typename Traits::template LocalView<ThisType, false>;
413 using SubControlVolume = typename Traits::SubControlVolume;
415 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
419 using DofMapper = typename Traits::ElementMapper;
420
423 static constexpr DiscretizationMethod discMethod{};
424
426 static constexpr int maxElementStencilSize = LocalView::maxNumElementScvfs*Traits::maxNumScvfNeighbors + 1;
427
429 using GridView = GV;
430
433 : ParentType(gridView)
434 {
435 // Check if the overlap size is what we expect
437 DUNE_THROW(Dune::InvalidStateException, "The cctpfa discretization method needs at least an overlap of 1 for parallel computations. "
438 << " Set the parameter \"Grid.Overlap\" in the input file.");
439
440 update_();
441 }
442
445 const DofMapper& dofMapper() const
446 { return this->elementMapper(); }
447
449 std::size_t numScv() const
450 {
451 return numScvs_;
452 }
453
455 std::size_t numScvf() const
456 {
457 return numScvf_;
458 }
459
461 std::size_t numBoundaryScvf() const
462 {
463 return numBoundaryScvf_;
464 }
465
467 std::size_t numDofs() const
468 { return this->gridView().size(0); }
469
471 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
472 void update()
473 {
474 ParentType::update();
475 update_();
476 }
477
479 void update(const GridView& gridView)
480 {
481 ParentType::update(gridView);
482 update_();
483 }
484
486 void update(GridView&& gridView)
487 {
488 ParentType::update(std::move(gridView));
489 update_();
490 }
491
492 const std::vector<GridIndexType>& scvfIndicesOfScv(GridIndexType scvIdx) const
493 { return scvfIndicesOfScv_[scvIdx]; }
494
496 const std::vector<NeighborVolVarIndices>& neighborVolVarIndices(GridIndexType scvIdx) const
497 { return neighborVolVarIndices_[scvIdx]; }
498
503 const ConnectivityMap &connectivityMap() const
504 { return connectivityMap_; }
505
506private:
507
508 void update_()
509 {
510 // clear local data
511 scvfIndicesOfScv_.clear();
512 neighborVolVarIndices_.clear();
513
514 // reserve memory or resize the containers
515 numScvs_ = numDofs();
516 numScvf_ = 0;
517 numBoundaryScvf_ = 0;
518 scvfIndicesOfScv_.resize(numScvs_);
519 neighborVolVarIndices_.resize(numScvs_);
520
521 // Build the SCV and SCV face
522 for (const auto& element : elements(this->gridView()))
523 {
524 const auto eIdx = this->elementMapper().index(element);
525
526 // the element-wise index sets for finite volume geometry
527 auto numLocalFaces = element.subEntities(1);
528 std::vector<GridIndexType> scvfsIndexSet;
529 std::vector<NeighborVolVarIndices> neighborVolVarIndexSet;
530 scvfsIndexSet.reserve(numLocalFaces);
531 neighborVolVarIndexSet.reserve(numLocalFaces);
532
533 // for network grids there might be multiple intersection with the same geometryInInside
534 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
535 std::vector<NeighborVolVarIndices> outsideIndices;
536 if (dim < dimWorld)
537 {
538 outsideIndices.resize(numLocalFaces);
539 for (const auto& intersection : intersections(this->gridView(), element))
540 {
541 if (intersection.neighbor())
542 {
543 const auto nIdx = this->elementMapper().index(intersection.outside());
544 outsideIndices[intersection.indexInInside()].push_back(nIdx);
545 }
546 }
547 }
548
549 for (const auto& intersection : intersections(this->gridView(), element))
550 {
551 // inner sub control volume faces (includes periodic boundaries)
552 if (intersection.neighbor())
553 {
554 // update the grid geometry if we have periodic boundaries
555 if (intersection.boundary())
556 this->setPeriodic();
557
558 if (dim == dimWorld)
559 {
560 scvfsIndexSet.push_back(numScvf_++);
561 const auto nIdx = this->elementMapper().index(intersection.outside());
562 neighborVolVarIndexSet.emplace_back(NeighborVolVarIndices({nIdx}));
563 }
564 // this is for network grids
565 // (will be optimized away of dim == dimWorld)
566 else
567 {
568 auto indexInInside = intersection.indexInInside();
569 // check if we already handled this facet
570 if (outsideIndices[indexInInside].empty())
571 continue;
572 else
573 {
574 scvfsIndexSet.push_back(numScvf_++);
575 neighborVolVarIndexSet.emplace_back(std::move(outsideIndices[indexInInside]));
576 outsideIndices[indexInInside].clear();
577 }
578 }
579 }
580 // boundary sub control volume faces
581 else if (intersection.boundary())
582 {
583 scvfsIndexSet.push_back(numScvf_++);
584 neighborVolVarIndexSet.emplace_back(NeighborVolVarIndices({static_cast<GridIndexType>(numScvs_ + numBoundaryScvf_++)}));
585 }
586 }
587
588 // store the sets of indices in the data container
589 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
590 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
591 }
592
593 // build the connectivity map for an efficient assembly
594 connectivityMap_.update(*this);
595 }
596
598 std::size_t numScvs_;
599 std::size_t numScvf_;
600 std::size_t numBoundaryScvf_;
601
603 ConnectivityMap connectivityMap_;
604
606 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
607 std::vector<std::vector<NeighborVolVarIndices>> neighborVolVarIndices_;
608};
609
610} // end namespace Dumux
611
612#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
Helper classes to compute the integration elements.
Base class for grid geometries.
Check the overlap size for different discretization methods.
The available discretization methods in Dumux.
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:51
GV GridView
export the grid view type
Definition: basegridgeometry.hh:66
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:52
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:55
static constexpr int maxNumScvfNeighbors
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:69
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:81
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:92
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:104
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:214
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:165
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:172
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:142
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:110
CCTpfaFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:125
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:106
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:108
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:160
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:205
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:148
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:192
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:179
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:218
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:154
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:199
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:112
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:186
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:138
The finite volume geometry (scvs and scvfs) for cell-centered TPFA models on a grid view This builds ...
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:393
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:411
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:492
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:419
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:445
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:415
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:486
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:413
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:455
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:503
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:472
CCTpfaFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:432
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:417
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:479
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:467
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:461
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:496
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/tpfa/fvgridgeometry.hh:449
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
Definition: method.hh:49
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...