3.5-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 <optional>
30#include <algorithm>
31#include <array>
32#include <vector>
33#include <utility>
34
35#include <dune/common/exceptions.hh>
37#include <dune/common/iteratorrange.hh>
39
40namespace Dumux {
41
51template<class GG, bool enableGridGeometryCache>
53
60template<class GG>
62{
64 using GridView = typename GG::GridView;
65 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
66
67public:
69 using Element = typename GridView::template Codim<0>::Entity;
71 using SubControlVolume = typename GG::SubControlVolume;
73 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
75 using GridGeometry = GG;
76
78 static constexpr std::size_t maxNumElementScvs = 1;
80 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
81
84 : gridGeometryPtr_(&gridGeometry) {}
85
88 const SubControlVolume& scv(GridIndexType scvIdx) const
89 {
90 return gridGeometry().scv(scvIdx);
91 }
92
95 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
96 {
97 return gridGeometry().scvf(scvfIdx);
98 }
99
102 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
103 {
104 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
105 }
106
112 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
113 scvs(const CCTpfaFVElementGeometry& fvGeometry)
114 {
116 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
117 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
118 }
119
125 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
127 {
128 const auto& g = fvGeometry.gridGeometry();
129 const auto scvIdx = fvGeometry.scvIndices_[0];
131 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
132 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
133 }
134
136 std::size_t numScv() const
137 {
138 return scvIndices_.size();
139 }
140
142 std::size_t numScvf() const
143 {
144 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
145 }
146
153 {
154 this->bindElement(element);
155 return std::move(*this);
156 }
157
158 void bind(const Element& element) &
159 {
160 this->bindElement(element);
161 }
162
169 {
170 this->bindElement(element);
171 return std::move(*this);
172 }
173
175 void bindElement(const Element& element) &
176 {
177 element_ = element;
178 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
179 }
180
182 bool isBound() const
183 { return static_cast<bool>(element_); }
184
186 const Element& element() const
187 { return *element_; }
188
191 { return *gridGeometryPtr_; }
192
194 bool hasBoundaryScvf() const
195 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
196
197private:
198
199 std::optional<Element> element_;
200 std::array<GridIndexType, 1> scvIndices_;
201 const GridGeometry* gridGeometryPtr_;
202};
203
209template<class GG>
211{
213 using GridView = typename GG::GridView;
214 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
215 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
216
217 static const int dim = GridView::dimension;
218 static const int dimWorld = GridView::dimensionworld;
219
220public:
222 using Element = typename GridView::template Codim<0>::Entity;
224 using SubControlVolume = typename GG::SubControlVolume;
226 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
228 using GridGeometry = GG;
230 static constexpr std::size_t maxNumElementScvs = 1;
232 static constexpr std::size_t maxNumElementScvfs = 2*dim;
233
236 : gridGeometryPtr_(&gridGeometry) {}
237
240 const SubControlVolume& scv(GridIndexType scvIdx) const
241 {
242 if (scvIdx == scvIndices_[0])
243 return scvs_[0];
244 else
245 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
246 }
247
250 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
251 {
252 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
253 if (it != scvfIndices_.end())
254 return scvfs_[std::distance(scvfIndices_.begin(), it)];
255 else
256 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
257 }
258
261 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
262 {
263 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
264 if (it != scvfIndices_.end())
265 {
266 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
267 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
268 }
269 else
270 {
271 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
272 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
273 if (localFlippedIndex < scvfs_.size())
274 return scvfs_[localFlippedIndex];
275 else
276 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
277 }
278 }
279
285 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
286 scvs(const ThisType& g)
287 {
288 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
289 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
290 }
291
297 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
298 scvfs(const ThisType& g)
299 {
300 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
301 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
302 }
303
305 std::size_t numScv() const
306 { return scvs_.size(); }
307
309 std::size_t numScvf() const
310 { return scvfs_.size(); }
311
318 {
319 this->bind_(element);
320 return std::move(*this);
321 }
322
323 void bind(const Element& element) &
324 {
325 this->bind_(element);
326 }
327
334 {
335 this->bindElement_(element);
336 return std::move(*this);
337 }
338
339 void bindElement(const Element& element) &
340 {
341 this->bindElement_(element);
342 }
343
345 bool isBound() const
346 { return static_cast<bool>(element_); }
347
349 const Element& element() const
350 { return *element_; }
351
354 { return *gridGeometryPtr_; }
355
357 bool hasBoundaryScvf() const
358 { return hasBoundaryScvf_; }
359
360private:
363 void bind_(const Element& element)
364 {
365 bindElement_(element);
366
367 neighborScvs_.reserve(element.subEntities(1));
368 neighborScvfIndices_.reserve(element.subEntities(1));
369 neighborScvfs_.reserve(element.subEntities(1));
370
371 std::vector<GridIndexType> handledNeighbors;
372 handledNeighbors.reserve(element.subEntities(1));
373
374 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
375 {
376 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
377 if (intersection.neighbor())
378 {
379 const auto outside = intersection.outside();
380 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
381
382 // make outside geometries only if not done yet (could happen on non-conforming grids)
383 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
384 {
385 makeNeighborGeometries(outside, outsideIdx);
386 handledNeighbors.push_back(outsideIdx);
387 }
388 }
389 }
390
391 // build flip index set for network, surface, and periodic grids
392 if (dim < dimWorld || gridGeometry().isPeriodic())
393 {
394 flippedScvfIndices_.resize(scvfs_.size());
395 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
396 {
397 const auto& scvf = scvfs_[localScvfIdx];
398 if (scvf.boundary())
399 continue;
400
401 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
402 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
403 {
404 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
405 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
406 {
407 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
408 {
409 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
410 break;
411 }
412 }
413 }
414 }
415
416 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
417 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
418 {
419 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
420 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
421 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
422 {
423 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
424 }
425 }
426 }
427 }
428
430 void bindElement_(const Element& element)
431 {
432 clear();
433 element_ = element;
434 scvfs_.reserve(element.subEntities(1));
435 scvfIndices_.reserve(element.subEntities(1));
436 makeElementGeometries(element);
437 }
438
439 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
440 {
441 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
442 {
443 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
444 {
445 return scvfs_.size() + localNeighborScvfIdx;
446 }
447 }
448
449 // go over all potential scvfs of the outside scv
450 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
451 {
452 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
453 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
454 {
455 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
456 {
457 return localOutsideScvfIdx;
458 }
459 }
460 }
461
462 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
463 }
464
466 void makeElementGeometries(const Element& element)
467 {
468 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
469
470 const auto eIdx = gridGeometry().elementMapper().index(element);
471 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
472 scvIndices_[0] = eIdx;
473
474 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
475 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
476
477 // for network grids there might be multiple intersection with the same geometryInInside
478 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
479 // here we keep track of them
480 std::vector<bool> handledScvf;
481 if (dim < dimWorld)
482 handledScvf.resize(element.subEntities(1), false);
483
484 int scvfCounter = 0;
485 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
486 {
487 if (dim < dimWorld)
488 if (handledScvf[intersection.indexInInside()])
489 continue;
490
491 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
492 if (intersection.neighbor() || intersection.boundary())
493 {
494 ScvfGridIndexStorage scvIndices;
495 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
496 scvIndices[0] = eIdx;
497 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
498
499 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
500 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
501
502 scvfs_.emplace_back(intersection,
503 intersection.geometry(),
504 scvFaceIndices[scvfCounter],
505 scvIndices,
506 onBoundary);
507 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
508 scvfCounter++;
509
510 // for surface and network grids mark that we handled this face
511 if (dim < dimWorld)
512 handledScvf[intersection.indexInInside()] = true;
513 }
514 }
515 }
516
518 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
519 {
520 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
521
522 // create the neighbor scv
523 neighborScvs_.emplace_back(element.geometry(), eIdx);
524 neighborScvIndices_.push_back(eIdx);
525
526 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
527 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
528
529 // for network grids there might be multiple intersection with the same geometryInInside
530 // we indentify those by the indexInInside for now (assumes conforming grids at branching facets)
531 // here we keep track of them
532 std::vector<bool> handledScvf;
533 if (dim < dimWorld)
534 handledScvf.resize(element.subEntities(1), false);
535
536 int scvfCounter = 0;
537 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
538 {
539 if (dim < dimWorld)
540 if (handledScvf[intersection.indexInInside()])
541 continue;
542
543 if (intersection.neighbor())
544 {
545 // this catches inner and periodic scvfs
546 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
547 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
548 {
549 // only create subcontrol faces where the outside element is the bound element
550 if (dim == dimWorld)
551 {
552 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
553 {
554 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
555 neighborScvfs_.emplace_back(intersection,
556 intersection.geometry(),
557 scvFaceIndices[scvfCounter],
558 scvIndices,
559 false);
560
561 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
562 }
563 }
564 // for network grids we can't use the intersection.outside() index as we can't assure that the
565 // first intersection with this indexInInside is the one that has our bound element as outside
566 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
567 // (will be optimized away for dim == dimWorld)
568 else
569 {
570 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
571 {
572 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
573 {
574 ScvfGridIndexStorage scvIndices;
575 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
576 scvIndices[0] = eIdx;
577 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
578 neighborScvfs_.emplace_back(intersection,
579 intersection.geometry(),
580 scvFaceIndices[scvfCounter],
581 scvIndices,
582 false);
583
584 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
585 break;
586 }
587 }
588 }
589
590 // for surface and network grids mark that we handled this face
591 if (dim < dimWorld)
592 handledScvf[intersection.indexInInside()] = true;
593 scvfCounter++;
594 }
595 }
596
597 // only increase counter for boundary intersections
598 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
599 else if (intersection.boundary())
600 {
601 if (dim < dimWorld)
602 handledScvf[intersection.indexInInside()] = true;
603 scvfCounter++;
604 }
605 }
606 }
607
608 const LocalIndexType findLocalIndex(const GridIndexType idx,
609 const std::vector<GridIndexType>& indices) const
610 {
611 auto it = std::find(indices.begin(), indices.end(), idx);
612 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
613 return std::distance(indices.begin(), it);
614
615 }
616
618 void clear()
619 {
620 scvfIndices_.clear();
621 scvfs_.clear();
622 flippedScvfIndices_.clear();
623
624 neighborScvIndices_.clear();
625 neighborScvfIndices_.clear();
626 neighborScvs_.clear();
627 neighborScvfs_.clear();
628 flippedNeighborScvfIndices_.clear();
629
630 hasBoundaryScvf_ = false;
631 }
632
633 std::optional<Element> element_;
634 const GridGeometry* gridGeometryPtr_;
635
636 // local storage after binding an element
637 std::array<GridIndexType, 1> scvIndices_;
638 std::array<SubControlVolume, 1> scvs_;
639
640 std::vector<GridIndexType> scvfIndices_;
641 std::vector<SubControlVolumeFace> scvfs_;
642 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
643
644 std::vector<GridIndexType> neighborScvIndices_;
645 std::vector<SubControlVolume> neighborScvs_;
646
647 std::vector<GridIndexType> neighborScvfIndices_;
648 std::vector<SubControlVolumeFace> neighborScvfs_;
649 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
650
651 bool hasBoundaryScvf_ = false;
652};
653
654} // end namespace Dumux
655
656#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.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:292
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:52
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:62
CCTpfaFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:168
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:69
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:186
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:158
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:83
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:126
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:136
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:194
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:75
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:190
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:142
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:182
CCTpfaFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:152
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:88
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:95
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:102
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:175
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:113
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:211
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:345
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:298
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:349
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:226
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:240
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:305
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:235
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:309
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:339
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:222
CCTpfaFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:317
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:353
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:323
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:261
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:250
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:224
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
CCTpfaFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:333
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:357
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82