3.6-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
197 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
198 {
199 assert(isBound());
200 return scv.geometry();
201 }
202
203 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
204 {
205 assert(isBound());
206 return scvf.geometry();
207 }
208
209private:
210
211 std::optional<Element> element_;
212 std::array<GridIndexType, 1> scvIndices_;
213 const GridGeometry* gridGeometryPtr_;
214};
215
221template<class GG>
223{
225 using GridView = typename GG::GridView;
226 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
227 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
228
229 static const int dim = GridView::dimension;
230 static const int dimWorld = GridView::dimensionworld;
231
232public:
234 using Element = typename GridView::template Codim<0>::Entity;
236 using SubControlVolume = typename GG::SubControlVolume;
238 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
240 using GridGeometry = GG;
242 static constexpr std::size_t maxNumElementScvs = 1;
244 static constexpr std::size_t maxNumElementScvfs = 2*dim;
245
248 : gridGeometryPtr_(&gridGeometry) {}
249
252 const SubControlVolume& scv(GridIndexType scvIdx) const
253 {
254 if (scvIdx == scvIndices_[0])
255 return scvs_[0];
256 else
257 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
258 }
259
262 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
263 {
264 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
265 if (it != scvfIndices_.end())
266 return scvfs_[std::distance(scvfIndices_.begin(), it)];
267 else
268 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
269 }
270
273 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
274 {
275 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
276 if (it != scvfIndices_.end())
277 {
278 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
279 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
280 }
281 else
282 {
283 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
284 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
285 if (localFlippedIndex < scvfs_.size())
286 return scvfs_[localFlippedIndex];
287 else
288 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
289 }
290 }
291
297 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
298 scvs(const ThisType& g)
299 {
300 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
301 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
302 }
303
309 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
310 scvfs(const ThisType& g)
311 {
312 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
313 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
314 }
315
317 std::size_t numScv() const
318 { return scvs_.size(); }
319
321 std::size_t numScvf() const
322 { return scvfs_.size(); }
323
330 {
331 this->bind_(element);
332 return std::move(*this);
333 }
334
335 void bind(const Element& element) &
336 {
337 this->bind_(element);
338 }
339
346 {
347 this->bindElement_(element);
348 return std::move(*this);
349 }
350
351 void bindElement(const Element& element) &
352 {
353 this->bindElement_(element);
354 }
355
357 bool isBound() const
358 { return static_cast<bool>(element_); }
359
361 const Element& element() const
362 { return *element_; }
363
366 { return *gridGeometryPtr_; }
367
369 bool hasBoundaryScvf() const
370 { return hasBoundaryScvf_; }
371
372 typename SubControlVolume::Traits::Geometry geometry (const SubControlVolume& scv) const
373 {
374 assert(isBound());
375 return scv.geometry();
376 }
377
378 typename SubControlVolumeFace::Traits::Geometry geometry (const SubControlVolumeFace& scvf) const
379 {
380 assert(isBound());
381 return scvf.geometry();
382 }
383
384private:
387 void bind_(const Element& element)
388 {
389 bindElement_(element);
390
391 neighborScvs_.reserve(element.subEntities(1));
392 neighborScvfIndices_.reserve(element.subEntities(1));
393 neighborScvfs_.reserve(element.subEntities(1));
394
395 std::vector<GridIndexType> handledNeighbors;
396 handledNeighbors.reserve(element.subEntities(1));
397
398 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
399 {
400 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
401 if (intersection.neighbor())
402 {
403 const auto outside = intersection.outside();
404 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
405
406 // make outside geometries only if not done yet (could happen on non-conforming grids)
407 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
408 {
409 makeNeighborGeometries(outside, outsideIdx);
410 handledNeighbors.push_back(outsideIdx);
411 }
412 }
413 }
414
415 // build flip index set for network, surface, and periodic grids
416 if (dim < dimWorld || gridGeometry().isPeriodic())
417 {
418 flippedScvfIndices_.resize(scvfs_.size());
419 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
420 {
421 const auto& scvf = scvfs_[localScvfIdx];
422 if (scvf.boundary())
423 continue;
424
425 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
426 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
427 {
428 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
429 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
430 {
431 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
432 {
433 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
434 break;
435 }
436 }
437 }
438 }
439
440 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
441 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
442 {
443 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
444 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
445 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
446 {
447 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
448 }
449 }
450 }
451 }
452
454 void bindElement_(const Element& element)
455 {
456 clear();
457 element_ = element;
458 scvfs_.reserve(element.subEntities(1));
459 scvfIndices_.reserve(element.subEntities(1));
460 makeElementGeometries(element);
461 }
462
463 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
464 {
465 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
466 {
467 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
468 {
469 return scvfs_.size() + localNeighborScvfIdx;
470 }
471 }
472
473 // go over all potential scvfs of the outside scv
474 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
475 {
476 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
477 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
478 {
479 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
480 {
481 return localOutsideScvfIdx;
482 }
483 }
484 }
485
486 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
487 }
488
490 void makeElementGeometries(const Element& element)
491 {
492 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
493
494 const auto eIdx = gridGeometry().elementMapper().index(element);
495 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
496 scvIndices_[0] = eIdx;
497
498 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
499 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
500
501 // for network grids there might be multiple intersection with the same geometryInInside
502 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
503 // here we keep track of them
504 std::vector<bool> handledScvf;
505 if (dim < dimWorld)
506 handledScvf.resize(element.subEntities(1), false);
507
508 int scvfCounter = 0;
509 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
510 {
511 if (dim < dimWorld)
512 if (handledScvf[intersection.indexInInside()])
513 continue;
514
515 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
516 if (intersection.neighbor() || intersection.boundary())
517 {
518 ScvfGridIndexStorage scvIndices;
519 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
520 scvIndices[0] = eIdx;
521 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
522
523 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
524 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
525
526 scvfs_.emplace_back(intersection,
527 intersection.geometry(),
528 scvFaceIndices[scvfCounter],
529 scvIndices,
530 onBoundary);
531 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
532 scvfCounter++;
533
534 // for surface and network grids mark that we handled this face
535 if (dim < dimWorld)
536 handledScvf[intersection.indexInInside()] = true;
537 }
538 }
539 }
540
542 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
543 {
544 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
545
546 // create the neighbor scv
547 neighborScvs_.emplace_back(element.geometry(), eIdx);
548 neighborScvIndices_.push_back(eIdx);
549
550 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
551 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
552
553 // for network grids there might be multiple intersection with the same geometryInInside
554 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
555 // here we keep track of them
556 std::vector<bool> handledScvf;
557 if (dim < dimWorld)
558 handledScvf.resize(element.subEntities(1), false);
559
560 int scvfCounter = 0;
561 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
562 {
563 if (dim < dimWorld)
564 if (handledScvf[intersection.indexInInside()])
565 continue;
566
567 if (intersection.neighbor())
568 {
569 // this catches inner and periodic scvfs
570 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
571 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
572 {
573 // only create subcontrol faces where the outside element is the bound element
574 if (dim == dimWorld)
575 {
576 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
577 {
578 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
579 neighborScvfs_.emplace_back(intersection,
580 intersection.geometry(),
581 scvFaceIndices[scvfCounter],
582 scvIndices,
583 false);
584
585 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
586 }
587 }
588 // for network grids we can't use the intersection.outside() index as we can't assure that the
589 // first intersection with this indexInInside is the one that has our bound element as outside
590 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
591 // (will be optimized away for dim == dimWorld)
592 else
593 {
594 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
595 {
596 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
597 {
598 ScvfGridIndexStorage scvIndices;
599 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
600 scvIndices[0] = eIdx;
601 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
602 neighborScvfs_.emplace_back(intersection,
603 intersection.geometry(),
604 scvFaceIndices[scvfCounter],
605 scvIndices,
606 false);
607
608 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
609 break;
610 }
611 }
612 }
613
614 // for surface and network grids mark that we handled this face
615 if (dim < dimWorld)
616 handledScvf[intersection.indexInInside()] = true;
617 scvfCounter++;
618 }
619 }
620
621 // only increase counter for boundary intersections
622 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
623 else if (intersection.boundary())
624 {
625 if (dim < dimWorld)
626 handledScvf[intersection.indexInInside()] = true;
627 scvfCounter++;
628 }
629 }
630 }
631
632 const LocalIndexType findLocalIndex(const GridIndexType idx,
633 const std::vector<GridIndexType>& indices) const
634 {
635 auto it = std::find(indices.begin(), indices.end(), idx);
636 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
637 return std::distance(indices.begin(), it);
638
639 }
640
642 void clear()
643 {
644 scvfIndices_.clear();
645 scvfs_.clear();
646 flippedScvfIndices_.clear();
647
648 neighborScvIndices_.clear();
649 neighborScvfIndices_.clear();
650 neighborScvs_.clear();
651 neighborScvfs_.clear();
652 flippedNeighborScvfIndices_.clear();
653
654 hasBoundaryScvf_ = false;
655 }
656
657 std::optional<Element> element_;
658 const GridGeometry* gridGeometryPtr_;
659
660 // local storage after binding an element
661 std::array<GridIndexType, 1> scvIndices_;
662 std::array<SubControlVolume, 1> scvs_;
663
664 std::vector<GridIndexType> scvfIndices_;
665 std::vector<SubControlVolumeFace> scvfs_;
666 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
667
668 std::vector<GridIndexType> neighborScvIndices_;
669 std::vector<SubControlVolume> neighborScvs_;
670
671 std::vector<GridIndexType> neighborScvfIndices_;
672 std::vector<SubControlVolumeFace> neighborScvfs_;
673 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
674
675 bool hasBoundaryScvf_ = false;
676};
677
678} // end namespace Dumux
679
680#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:294
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:197
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
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:203
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:223
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:357
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:310
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:361
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:238
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:252
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:317
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:247
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:321
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:372
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:351
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:234
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:329
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:378
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:365
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:335
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:273
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:262
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:236
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:298
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:345
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:240
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:369
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82