3.1-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 using FVGridGeometry [[deprecated("Use GridGeometry instead. Will be removed after 3.1!")]] = GridGeometry;
74
76 static constexpr std::size_t maxNumElementScvs = 1;
78 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
79
82 : gridGeometryPtr_(&gridGeometry) {}
83
86 const SubControlVolume& scv(GridIndexType scvIdx) const
87 {
88 return gridGeometry().scv(scvIdx);
89 }
90
93 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
94 {
95 return gridGeometry().scvf(scvfIdx);
96 }
97
100 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
101 {
102 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
103 }
104
110 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
111 scvs(const CCTpfaFVElementGeometry& fvGeometry)
112 {
114 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
115 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
116 }
117
123 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
125 {
126 const auto& g = fvGeometry.gridGeometry();
127 const auto scvIdx = fvGeometry.scvIndices_[0];
129 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
130 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
131 }
132
134 std::size_t numScv() const
135 {
136 return scvIndices_.size();
137 }
138
140 std::size_t numScvf() const
141 {
142 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
143 }
144
146 void bind(const Element& element)
147 {
148 this->bindElement(element);
149 }
150
152 void bindElement(const Element& element)
153 {
154 elementPtr_ = &element;
155 scvIndices_[0] = gridGeometry().elementMapper().index(*elementPtr_);
156 }
157
160 { return *gridGeometryPtr_; }
161 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
163 { return gridGeometry(); }
164
166 bool hasBoundaryScvf() const
167 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
168
169private:
170
171 const Element* elementPtr_;
172 std::array<GridIndexType, 1> scvIndices_;
173 const GridGeometry* gridGeometryPtr_;
174};
175
181template<class GG>
183{
185 using GridView = typename GG::GridView;
186 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
187 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
188 using Element = typename GridView::template Codim<0>::Entity;
189
190 static const int dim = GridView::dimension;
191 static const int dimWorld = GridView::dimensionworld;
192
193public:
195 using SubControlVolume = typename GG::SubControlVolume;
197 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
199 using GridGeometry = GG;
200 using FVGridGeometry [[deprecated("Use GridGeometry instead. Will be removed after 3.1!")]] = GridGeometry;
202 static constexpr std::size_t maxNumElementScvs = 1;
204 static constexpr std::size_t maxNumElementScvfs = 2*dim;
205
208 : gridGeometryPtr_(&gridGeometry) {}
209
212 const SubControlVolume& scv(GridIndexType scvIdx) const
213 {
214 if (scvIdx == scvIndices_[0])
215 return scvs_[0];
216 else
217 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
218 }
219
222 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
223 {
224 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
225 if (it != scvfIndices_.end())
226 return scvfs_[std::distance(scvfIndices_.begin(), it)];
227 else
228 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
229 }
230
233 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
234 {
235 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
236 if (it != scvfIndices_.end())
237 {
238 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
239 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
240 }
241 else
242 {
243 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
244 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
245 if (localFlippedIndex < scvfs_.size())
246 return scvfs_[localFlippedIndex];
247 else
248 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
249 }
250 }
251
257 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
258 scvs(const ThisType& g)
259 {
260 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
261 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
262 }
263
269 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
270 scvfs(const ThisType& g)
271 {
272 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
273 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
274 }
275
277 std::size_t numScv() const
278 { return scvs_.size(); }
279
281 std::size_t numScvf() const
282 { return scvfs_.size(); }
283
286 void bind(const Element& element)
287 {
288 bindElement(element);
289
290 neighborScvs_.reserve(element.subEntities(1));
291 neighborScvfIndices_.reserve(element.subEntities(1));
292 neighborScvfs_.reserve(element.subEntities(1));
293
294 std::vector<GridIndexType> handledNeighbors;
295 handledNeighbors.reserve(element.subEntities(1));
296
297 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
298 {
299 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
300 if (intersection.neighbor())
301 {
302 const auto outside = intersection.outside();
303 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
304
305 // make outside geometries only if not done yet (could happen on non-conforming grids)
306 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
307 {
308 makeNeighborGeometries(outside, outsideIdx);
309 handledNeighbors.push_back(outsideIdx);
310 }
311 }
312 }
313
314 // build flip index set for network, surface, and periodic grids
315 if (dim < dimWorld || gridGeometry().isPeriodic())
316 {
317 flippedScvfIndices_.resize(scvfs_.size());
318 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
319 {
320 const auto& scvf = scvfs_[localScvfIdx];
321 if (scvf.boundary())
322 continue;
323
324 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
325 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
326 {
327 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
328 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
329 {
330 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
331 {
332 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
333 break;
334 }
335 }
336 }
337 }
338
339 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
340 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
341 {
342 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
343 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
344 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
345 {
346 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
347 }
348 }
349 }
350 }
351
353 void bindElement(const Element& element)
354 {
355 clear();
356 elementPtr_ = &element;
357 scvfs_.reserve(element.subEntities(1));
358 scvfIndices_.reserve(element.subEntities(1));
359 makeElementGeometries(element);
360 }
361
364 { return *gridGeometryPtr_; }
365 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
367 { return gridGeometry(); }
368
370 bool hasBoundaryScvf() const
371 { return hasBoundaryScvf_; }
372
373private:
374
375 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
376 {
377 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
378 {
379 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
380 {
381 return scvfs_.size() + localNeighborScvfIdx;
382 }
383 }
384
385 // go over all potential scvfs of the outside scv
386 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
387 {
388 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
389 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
390 {
391 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
392 {
393 return localOutsideScvfIdx;
394 }
395 }
396 }
397
398 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
399 }
400
402 void makeElementGeometries(const Element& element)
403 {
404 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
405
406 const auto eIdx = gridGeometry().elementMapper().index(element);
407 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
408 scvIndices_[0] = eIdx;
409
410 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
411 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
412
413 // for network grids there might be multiple intersection with the same geometryInInside
414 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
415 // here we keep track of them
416 std::vector<bool> handledScvf;
417 if (dim < dimWorld)
418 handledScvf.resize(element.subEntities(1), false);
419
420 int scvfCounter = 0;
421 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
422 {
423 if (dim < dimWorld)
424 if (handledScvf[intersection.indexInInside()])
425 continue;
426
427 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
428 if (intersection.neighbor() || intersection.boundary())
429 {
430 ScvfGridIndexStorage scvIndices;
431 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
432 scvIndices[0] = eIdx;
433 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
434
435 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
436 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
437
438 scvfs_.emplace_back(intersection,
439 intersection.geometry(),
440 scvFaceIndices[scvfCounter],
441 scvIndices,
442 onBoundary);
443 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
444 scvfCounter++;
445
446 // for surface and network grids mark that we handled this face
447 if (dim < dimWorld)
448 handledScvf[intersection.indexInInside()] = true;
449 }
450 }
451 }
452
454 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
455 {
456 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
457
458 // create the neighbor scv
459 neighborScvs_.emplace_back(element.geometry(), eIdx);
460 neighborScvIndices_.push_back(eIdx);
461
462 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
463 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
464
465 // for network grids there might be multiple intersection with the same geometryInInside
466 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
467 // here we keep track of them
468 std::vector<bool> handledScvf;
469 if (dim < dimWorld)
470 handledScvf.resize(element.subEntities(1), false);
471
472 int scvfCounter = 0;
473 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
474 {
475 if (dim < dimWorld)
476 if (handledScvf[intersection.indexInInside()])
477 continue;
478
479 if (intersection.neighbor())
480 {
481 // this catches inner and periodic scvfs
482 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
483 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
484 {
485 // only create subcontrol faces where the outside element is the bound element
486 if (dim == dimWorld)
487 {
488 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*elementPtr_))
489 {
490 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
491 neighborScvfs_.emplace_back(intersection,
492 intersection.geometry(),
493 scvFaceIndices[scvfCounter],
494 scvIndices,
495 false);
496
497 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
498 }
499 }
500 // for network grids we can't use the intersection.outside() index as we can't assure that the
501 // first intersection with this indexInInside is the one that has our bound element as outside
502 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
503 // (will be optimized away for dim == dimWorld)
504 else
505 {
506 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
507 {
508 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*elementPtr_))
509 {
510 ScvfGridIndexStorage scvIndices;
511 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
512 scvIndices[0] = eIdx;
513 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
514 neighborScvfs_.emplace_back(intersection,
515 intersection.geometry(),
516 scvFaceIndices[scvfCounter],
517 scvIndices,
518 false);
519
520 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
521 break;
522 }
523 }
524 }
525
526 // for surface and network grids mark that we handled this face
527 if (dim < dimWorld)
528 handledScvf[intersection.indexInInside()] = true;
529 scvfCounter++;
530 }
531 }
532
533 // only increase counter for boundary intersections
534 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
535 else if (intersection.boundary())
536 {
537 if (dim < dimWorld)
538 handledScvf[intersection.indexInInside()] = true;
539 scvfCounter++;
540 }
541 }
542 }
543
544 const LocalIndexType findLocalIndex(const GridIndexType idx,
545 const std::vector<GridIndexType>& indices) const
546 {
547 auto it = std::find(indices.begin(), indices.end(), idx);
548 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
549 return std::distance(indices.begin(), it);
550
551 }
552
554 void clear()
555 {
556 scvfIndices_.clear();
557 scvfs_.clear();
558 flippedScvfIndices_.clear();
559
560 neighborScvIndices_.clear();
561 neighborScvfIndices_.clear();
562 neighborScvs_.clear();
563 neighborScvfs_.clear();
564 flippedNeighborScvfIndices_.clear();
565
566 hasBoundaryScvf_ = false;
567 }
568
569 const Element* elementPtr_;
570 const GridGeometry* gridGeometryPtr_;
571
572 // local storage after binding an element
573 std::array<GridIndexType, 1> scvIndices_;
574 std::array<SubControlVolume, 1> scvs_;
575
576 std::vector<GridIndexType> scvfIndices_;
577 std::vector<SubControlVolumeFace> scvfs_;
578 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
579
580 std::vector<GridIndexType> neighborScvIndices_;
581 std::vector<SubControlVolume> neighborScvs_;
582
583 std::vector<GridIndexType> neighborScvfIndices_;
584 std::vector<SubControlVolumeFace> neighborScvfs_;
585 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
586
587 bool hasBoundaryScvf_ = false;
588};
589
590} // end namespace Dumux
591
592#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.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
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:81
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:124
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:134
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:166
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:146
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:159
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:140
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:86
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:152
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:93
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:100
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:111
const GridGeometry & fvGridGeometry() const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:162
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:183
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:270
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:197
const GridGeometry & fvGridGeometry() const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:366
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:212
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:277
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:207
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:353
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:281
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:363
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:233
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:222
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:200
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:195
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:258
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:199
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:370
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82