Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 <optional>
30#include <algorithm>
31#include <array>
32#include <vector>
33
34#include <dune/common/exceptions.hh>
36#include <dune/common/iteratorrange.hh>
38
39namespace Dumux {
40
50template<class GG, bool enableGridGeometryCache>
52
59template<class GG>
61{
63 using GridView = typename GG::GridView;
64 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
65
66public:
68 using Element = typename GridView::template Codim<0>::Entity;
70 using SubControlVolume = typename GG::SubControlVolume;
72 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
74 using GridGeometry = GG;
75
77 static constexpr std::size_t maxNumElementScvs = 1;
79 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
80
83 : gridGeometryPtr_(&gridGeometry) {}
84
87 const SubControlVolume& scv(GridIndexType scvIdx) const
88 {
89 return gridGeometry().scv(scvIdx);
90 }
91
94 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
95 {
96 return gridGeometry().scvf(scvfIdx);
97 }
98
101 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
102 {
103 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
104 }
105
111 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
112 scvs(const CCTpfaFVElementGeometry& fvGeometry)
113 {
115 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
116 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
117 }
118
124 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
126 {
127 const auto& g = fvGeometry.gridGeometry();
128 const auto scvIdx = fvGeometry.scvIndices_[0];
130 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
131 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
132 }
133
135 std::size_t numScv() const
136 {
137 return scvIndices_.size();
138 }
139
141 std::size_t numScvf() const
142 {
143 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
144 }
145
147 void bind(const Element& element)
148 {
149 this->bindElement(element);
150 }
151
153 void bindElement(const Element& element)
154 {
155 element_ = element;
156 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
157 }
158
160 bool isBound() const
161 { return static_cast<bool>(element_); }
162
164 const Element& element() const
165 { return *element_; }
166
169 { return *gridGeometryPtr_; }
170
172 bool hasBoundaryScvf() const
173 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
174
175private:
176
177 std::optional<Element> element_;
178 std::array<GridIndexType, 1> scvIndices_;
179 const GridGeometry* gridGeometryPtr_;
180};
181
187template<class GG>
189{
191 using GridView = typename GG::GridView;
192 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
193 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
194
195 static const int dim = GridView::dimension;
196 static const int dimWorld = GridView::dimensionworld;
197
198public:
200 using Element = typename GridView::template Codim<0>::Entity;
202 using SubControlVolume = typename GG::SubControlVolume;
204 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
206 using GridGeometry = GG;
208 static constexpr std::size_t maxNumElementScvs = 1;
210 static constexpr std::size_t maxNumElementScvfs = 2*dim;
211
214 : gridGeometryPtr_(&gridGeometry) {}
215
218 const SubControlVolume& scv(GridIndexType scvIdx) const
219 {
220 if (scvIdx == scvIndices_[0])
221 return scvs_[0];
222 else
223 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
224 }
225
228 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
229 {
230 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
231 if (it != scvfIndices_.end())
232 return scvfs_[std::distance(scvfIndices_.begin(), it)];
233 else
234 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
235 }
236
239 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
240 {
241 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
242 if (it != scvfIndices_.end())
243 {
244 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
245 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
246 }
247 else
248 {
249 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
250 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
251 if (localFlippedIndex < scvfs_.size())
252 return scvfs_[localFlippedIndex];
253 else
254 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
255 }
256 }
257
263 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
264 scvs(const ThisType& g)
265 {
266 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
267 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
268 }
269
275 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
276 scvfs(const ThisType& g)
277 {
278 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
279 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
280 }
281
283 std::size_t numScv() const
284 { return scvs_.size(); }
285
287 std::size_t numScvf() const
288 { return scvfs_.size(); }
289
292 void bind(const Element& element)
293 {
294 bindElement(element);
295
296 neighborScvs_.reserve(element.subEntities(1));
297 neighborScvfIndices_.reserve(element.subEntities(1));
298 neighborScvfs_.reserve(element.subEntities(1));
299
300 std::vector<GridIndexType> handledNeighbors;
301 handledNeighbors.reserve(element.subEntities(1));
302
303 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
304 {
305 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
306 if (intersection.neighbor())
307 {
308 const auto outside = intersection.outside();
309 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
310
311 // make outside geometries only if not done yet (could happen on non-conforming grids)
312 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
313 {
314 makeNeighborGeometries(outside, outsideIdx);
315 handledNeighbors.push_back(outsideIdx);
316 }
317 }
318 }
319
320 // build flip index set for network, surface, and periodic grids
321 if (dim < dimWorld || gridGeometry().isPeriodic())
322 {
323 flippedScvfIndices_.resize(scvfs_.size());
324 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
325 {
326 const auto& scvf = scvfs_[localScvfIdx];
327 if (scvf.boundary())
328 continue;
329
330 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
331 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
332 {
333 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
334 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
335 {
336 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
337 {
338 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
339 break;
340 }
341 }
342 }
343 }
344
345 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
346 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
347 {
348 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
349 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
350 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
351 {
352 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
353 }
354 }
355 }
356 }
357
359 void bindElement(const Element& element)
360 {
361 clear();
362 element_ = element;
363 scvfs_.reserve(element.subEntities(1));
364 scvfIndices_.reserve(element.subEntities(1));
365 makeElementGeometries(element);
366 }
367
369 bool isBound() const
370 { return static_cast<bool>(element_); }
371
373 const Element& element() const
374 { return *element_; }
375
378 { return *gridGeometryPtr_; }
379
381 bool hasBoundaryScvf() const
382 { return hasBoundaryScvf_; }
383
384private:
385
386 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
387 {
388 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
389 {
390 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
391 {
392 return scvfs_.size() + localNeighborScvfIdx;
393 }
394 }
395
396 // go over all potential scvfs of the outside scv
397 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
398 {
399 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
400 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
401 {
402 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
403 {
404 return localOutsideScvfIdx;
405 }
406 }
407 }
408
409 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
410 }
411
413 void makeElementGeometries(const Element& element)
414 {
415 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
416
417 const auto eIdx = gridGeometry().elementMapper().index(element);
418 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
419 scvIndices_[0] = eIdx;
420
421 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
422 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
423
424 // for network grids there might be multiple intersection with the same geometryInInside
425 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
426 // here we keep track of them
427 std::vector<bool> handledScvf;
428 if (dim < dimWorld)
429 handledScvf.resize(element.subEntities(1), false);
430
431 int scvfCounter = 0;
432 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
433 {
434 if (dim < dimWorld)
435 if (handledScvf[intersection.indexInInside()])
436 continue;
437
438 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
439 if (intersection.neighbor() || intersection.boundary())
440 {
441 ScvfGridIndexStorage scvIndices;
442 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
443 scvIndices[0] = eIdx;
444 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
445
446 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
447 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
448
449 scvfs_.emplace_back(intersection,
450 intersection.geometry(),
451 scvFaceIndices[scvfCounter],
452 scvIndices,
453 onBoundary);
454 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
455 scvfCounter++;
456
457 // for surface and network grids mark that we handled this face
458 if (dim < dimWorld)
459 handledScvf[intersection.indexInInside()] = true;
460 }
461 }
462 }
463
465 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
466 {
467 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
468
469 // create the neighbor scv
470 neighborScvs_.emplace_back(element.geometry(), eIdx);
471 neighborScvIndices_.push_back(eIdx);
472
473 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
474 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
475
476 // for network grids there might be multiple intersection with the same geometryInInside
477 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
478 // here we keep track of them
479 std::vector<bool> handledScvf;
480 if (dim < dimWorld)
481 handledScvf.resize(element.subEntities(1), false);
482
483 int scvfCounter = 0;
484 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
485 {
486 if (dim < dimWorld)
487 if (handledScvf[intersection.indexInInside()])
488 continue;
489
490 if (intersection.neighbor())
491 {
492 // this catches inner and periodic scvfs
493 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
494 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
495 {
496 // only create subcontrol faces where the outside element is the bound element
497 if (dim == dimWorld)
498 {
499 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
500 {
501 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
502 neighborScvfs_.emplace_back(intersection,
503 intersection.geometry(),
504 scvFaceIndices[scvfCounter],
505 scvIndices,
506 false);
507
508 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
509 }
510 }
511 // for network grids we can't use the intersection.outside() index as we can't assure that the
512 // first intersection with this indexInInside is the one that has our bound element as outside
513 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
514 // (will be optimized away for dim == dimWorld)
515 else
516 {
517 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
518 {
519 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
520 {
521 ScvfGridIndexStorage scvIndices;
522 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
523 scvIndices[0] = eIdx;
524 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
525 neighborScvfs_.emplace_back(intersection,
526 intersection.geometry(),
527 scvFaceIndices[scvfCounter],
528 scvIndices,
529 false);
530
531 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
532 break;
533 }
534 }
535 }
536
537 // for surface and network grids mark that we handled this face
538 if (dim < dimWorld)
539 handledScvf[intersection.indexInInside()] = true;
540 scvfCounter++;
541 }
542 }
543
544 // only increase counter for boundary intersections
545 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
546 else if (intersection.boundary())
547 {
548 if (dim < dimWorld)
549 handledScvf[intersection.indexInInside()] = true;
550 scvfCounter++;
551 }
552 }
553 }
554
555 const LocalIndexType findLocalIndex(const GridIndexType idx,
556 const std::vector<GridIndexType>& indices) const
557 {
558 auto it = std::find(indices.begin(), indices.end(), idx);
559 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
560 return std::distance(indices.begin(), it);
561
562 }
563
565 void clear()
566 {
567 scvfIndices_.clear();
568 scvfs_.clear();
569 flippedScvfIndices_.clear();
570
571 neighborScvIndices_.clear();
572 neighborScvfIndices_.clear();
573 neighborScvs_.clear();
574 neighborScvfs_.clear();
575 flippedNeighborScvfIndices_.clear();
576
577 hasBoundaryScvf_ = false;
578 }
579
580 std::optional<Element> element_;
581 const GridGeometry* gridGeometryPtr_;
582
583 // local storage after binding an element
584 std::array<GridIndexType, 1> scvIndices_;
585 std::array<SubControlVolume, 1> scvs_;
586
587 std::vector<GridIndexType> scvfIndices_;
588 std::vector<SubControlVolumeFace> scvfs_;
589 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
590
591 std::vector<GridIndexType> neighborScvIndices_;
592 std::vector<SubControlVolume> neighborScvs_;
593
594 std::vector<GridIndexType> neighborScvfIndices_;
595 std::vector<SubControlVolumeFace> neighborScvfs_;
596 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
597
598 bool hasBoundaryScvf_ = false;
599};
600
601} // end namespace Dumux
602
603#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.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
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:51
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:61
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:72
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:164
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:82
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:125
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:135
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:172
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:147
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:74
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:168
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:141
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:160
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:87
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:153
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:94
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:101
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:70
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:112
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:189
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:369
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:276
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:373
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:204
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:218
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:283
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:213
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:359
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:287
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:200
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:377
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:239
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:292
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:264
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:206
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:381
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82