3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/cellcentered/tpfa/fvelementgeometry.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_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
28
29#include <algorithm>
30#include <array>
31#include <vector>
32
33#include <dune/common/exceptions.hh>
35#include <dune/common/iteratorrange.hh>
37
38namespace Dumux {
39
49template<class GG, bool enableGridGeometryCache>
51
58template<class GG>
60{
62 using GridView = typename GG::GridView;
63 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
64 using Element = typename GridView::template Codim<0>::Entity;
65
66public:
68 using SubControlVolume = typename GG::SubControlVolume;
70 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
72 using GridGeometry = GG;
73
75 static constexpr std::size_t maxNumElementScvs = 1;
77 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
78
81 : gridGeometryPtr_(&gridGeometry) {}
82
85 const SubControlVolume& scv(GridIndexType scvIdx) const
86 {
87 return gridGeometry().scv(scvIdx);
88 }
89
92 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
93 {
94 return gridGeometry().scvf(scvfIdx);
95 }
96
99 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
100 {
101 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
102 }
103
109 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
110 scvs(const CCTpfaFVElementGeometry& fvGeometry)
111 {
113 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
114 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
115 }
116
122 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
124 {
125 const auto& g = fvGeometry.gridGeometry();
126 const auto scvIdx = fvGeometry.scvIndices_[0];
128 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
129 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
130 }
131
133 std::size_t numScv() const
134 {
135 return scvIndices_.size();
136 }
137
139 std::size_t numScvf() const
140 {
141 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
142 }
143
145 void bind(const Element& element)
146 {
147 this->bindElement(element);
148 }
149
151 void bindElement(const Element& element)
152 {
153 elementPtr_ = &element;
154 scvIndices_[0] = gridGeometry().elementMapper().index(*elementPtr_);
155 }
156
159 { return *gridGeometryPtr_; }
160
162 bool hasBoundaryScvf() const
163 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
164
165private:
166
167 const Element* elementPtr_;
168 std::array<GridIndexType, 1> scvIndices_;
169 const GridGeometry* gridGeometryPtr_;
170};
171
177template<class GG>
179{
181 using GridView = typename GG::GridView;
182 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
183 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
184 using Element = typename GridView::template Codim<0>::Entity;
185
186 static const int dim = GridView::dimension;
187 static const int dimWorld = GridView::dimensionworld;
188
189public:
191 using SubControlVolume = typename GG::SubControlVolume;
193 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
195 using GridGeometry = GG;
197 static constexpr std::size_t maxNumElementScvs = 1;
199 static constexpr std::size_t maxNumElementScvfs = 2*dim;
200
203 : gridGeometryPtr_(&gridGeometry) {}
204
207 const SubControlVolume& scv(GridIndexType scvIdx) const
208 {
209 if (scvIdx == scvIndices_[0])
210 return scvs_[0];
211 else
212 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
213 }
214
217 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
218 {
219 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
220 if (it != scvfIndices_.end())
221 return scvfs_[std::distance(scvfIndices_.begin(), it)];
222 else
223 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
224 }
225
228 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
229 {
230 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
231 if (it != scvfIndices_.end())
232 {
233 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
234 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
235 }
236 else
237 {
238 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
239 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
240 if (localFlippedIndex < scvfs_.size())
241 return scvfs_[localFlippedIndex];
242 else
243 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
244 }
245 }
246
252 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
253 scvs(const ThisType& g)
254 {
255 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
256 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
257 }
258
264 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
265 scvfs(const ThisType& g)
266 {
267 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
268 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
269 }
270
272 std::size_t numScv() const
273 { return scvs_.size(); }
274
276 std::size_t numScvf() const
277 { return scvfs_.size(); }
278
281 void bind(const Element& element)
282 {
283 bindElement(element);
284
285 neighborScvs_.reserve(element.subEntities(1));
286 neighborScvfIndices_.reserve(element.subEntities(1));
287 neighborScvfs_.reserve(element.subEntities(1));
288
289 std::vector<GridIndexType> handledNeighbors;
290 handledNeighbors.reserve(element.subEntities(1));
291
292 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
293 {
294 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
295 if (intersection.neighbor())
296 {
297 const auto outside = intersection.outside();
298 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
299
300 // make outside geometries only if not done yet (could happen on non-conforming grids)
301 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
302 {
303 makeNeighborGeometries(outside, outsideIdx);
304 handledNeighbors.push_back(outsideIdx);
305 }
306 }
307 }
308
309 // build flip index set for network, surface, and periodic grids
310 if (dim < dimWorld || gridGeometry().isPeriodic())
311 {
312 flippedScvfIndices_.resize(scvfs_.size());
313 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
314 {
315 const auto& scvf = scvfs_[localScvfIdx];
316 if (scvf.boundary())
317 continue;
318
319 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
320 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
321 {
322 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
323 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
324 {
325 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
326 {
327 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
328 break;
329 }
330 }
331 }
332 }
333
334 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
335 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
336 {
337 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
338 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
339 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
340 {
341 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
342 }
343 }
344 }
345 }
346
348 void bindElement(const Element& element)
349 {
350 clear();
351 elementPtr_ = &element;
352 scvfs_.reserve(element.subEntities(1));
353 scvfIndices_.reserve(element.subEntities(1));
354 makeElementGeometries(element);
355 }
356
359 { return *gridGeometryPtr_; }
360
362 bool hasBoundaryScvf() const
363 { return hasBoundaryScvf_; }
364
365private:
366
367 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
368 {
369 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
370 {
371 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
372 {
373 return scvfs_.size() + localNeighborScvfIdx;
374 }
375 }
376
377 // go over all potential scvfs of the outside scv
378 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
379 {
380 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
381 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
382 {
383 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
384 {
385 return localOutsideScvfIdx;
386 }
387 }
388 }
389
390 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
391 }
392
394 void makeElementGeometries(const Element& element)
395 {
396 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
397
398 const auto eIdx = gridGeometry().elementMapper().index(element);
399 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
400 scvIndices_[0] = eIdx;
401
402 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
403 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
404
405 // for network grids there might be multiple intersection with the same geometryInInside
406 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
407 // here we keep track of them
408 std::vector<bool> handledScvf;
409 if (dim < dimWorld)
410 handledScvf.resize(element.subEntities(1), false);
411
412 int scvfCounter = 0;
413 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
414 {
415 if (dim < dimWorld)
416 if (handledScvf[intersection.indexInInside()])
417 continue;
418
419 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
420 if (intersection.neighbor() || intersection.boundary())
421 {
422 ScvfGridIndexStorage scvIndices;
423 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
424 scvIndices[0] = eIdx;
425 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
426
427 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
428 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
429
430 scvfs_.emplace_back(intersection,
431 intersection.geometry(),
432 scvFaceIndices[scvfCounter],
433 scvIndices,
434 onBoundary);
435 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
436 scvfCounter++;
437
438 // for surface and network grids mark that we handled this face
439 if (dim < dimWorld)
440 handledScvf[intersection.indexInInside()] = true;
441 }
442 }
443 }
444
446 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
447 {
448 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
449
450 // create the neighbor scv
451 neighborScvs_.emplace_back(element.geometry(), eIdx);
452 neighborScvIndices_.push_back(eIdx);
453
454 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
455 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
456
457 // for network grids there might be multiple intersection with the same geometryInInside
458 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
459 // here we keep track of them
460 std::vector<bool> handledScvf;
461 if (dim < dimWorld)
462 handledScvf.resize(element.subEntities(1), false);
463
464 int scvfCounter = 0;
465 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
466 {
467 if (dim < dimWorld)
468 if (handledScvf[intersection.indexInInside()])
469 continue;
470
471 if (intersection.neighbor())
472 {
473 // this catches inner and periodic scvfs
474 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
475 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
476 {
477 // only create subcontrol faces where the outside element is the bound element
478 if (dim == dimWorld)
479 {
480 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*elementPtr_))
481 {
482 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
483 neighborScvfs_.emplace_back(intersection,
484 intersection.geometry(),
485 scvFaceIndices[scvfCounter],
486 scvIndices,
487 false);
488
489 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
490 }
491 }
492 // for network grids we can't use the intersection.outside() index as we can't assure that the
493 // first intersection with this indexInInside is the one that has our bound element as outside
494 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
495 // (will be optimized away for dim == dimWorld)
496 else
497 {
498 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
499 {
500 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*elementPtr_))
501 {
502 ScvfGridIndexStorage scvIndices;
503 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
504 scvIndices[0] = eIdx;
505 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
506 neighborScvfs_.emplace_back(intersection,
507 intersection.geometry(),
508 scvFaceIndices[scvfCounter],
509 scvIndices,
510 false);
511
512 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
513 break;
514 }
515 }
516 }
517
518 // for surface and network grids mark that we handled this face
519 if (dim < dimWorld)
520 handledScvf[intersection.indexInInside()] = true;
521 scvfCounter++;
522 }
523 }
524
525 // only increase counter for boundary intersections
526 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
527 else if (intersection.boundary())
528 {
529 if (dim < dimWorld)
530 handledScvf[intersection.indexInInside()] = true;
531 scvfCounter++;
532 }
533 }
534 }
535
536 const LocalIndexType findLocalIndex(const GridIndexType idx,
537 const std::vector<GridIndexType>& indices) const
538 {
539 auto it = std::find(indices.begin(), indices.end(), idx);
540 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
541 return std::distance(indices.begin(), it);
542
543 }
544
546 void clear()
547 {
548 scvfIndices_.clear();
549 scvfs_.clear();
550 flippedScvfIndices_.clear();
551
552 neighborScvIndices_.clear();
553 neighborScvfIndices_.clear();
554 neighborScvs_.clear();
555 neighborScvfs_.clear();
556 flippedNeighborScvfIndices_.clear();
557
558 hasBoundaryScvf_ = false;
559 }
560
561 const Element* elementPtr_;
562 const GridGeometry* gridGeometryPtr_;
563
564 // local storage after binding an element
565 std::array<GridIndexType, 1> scvIndices_;
566 std::array<SubControlVolume, 1> scvs_;
567
568 std::vector<GridIndexType> scvfIndices_;
569 std::vector<SubControlVolumeFace> scvfs_;
570 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
571
572 std::vector<GridIndexType> neighborScvIndices_;
573 std::vector<SubControlVolume> neighborScvs_;
574
575 std::vector<GridIndexType> neighborScvfIndices_;
576 std::vector<SubControlVolumeFace> neighborScvfs_;
577 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
578
579 bool hasBoundaryScvf_ = false;
580};
581
582} // end namespace Dumux
583
584#endif
Defines the index types used for grid and local indices.
Class providing iterators over sub control volumes and sub control volume faces of an element.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:50
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:60
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:70
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:80
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:123
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:133
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:162
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:145
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:72
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:158
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:139
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:85
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:151
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:92
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:99
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:110
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:179
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:265
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:193
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:207
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:272
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:348
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:276
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:358
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:217
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:281
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:191
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:253
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:195
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:362
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82