version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <algorithm>
19#include <array>
20#include <vector>
21#include <utility>
22
23#include <dune/common/exceptions.hh>
25#include <dune/common/iteratorrange.hh>
27
28namespace Dumux {
29
30namespace Detail::Tpfa {
31
32template<class GridIndexType>
33auto findLocalIndex(const GridIndexType idx,
34 const std::vector<GridIndexType>& indices)
35{
36 auto it = std::find(indices.begin(), indices.end(), idx);
37 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
38 return std::distance(indices.begin(), it);
39}
40
41} // end namespace Detail::Tpfa
42
52template<class GG, bool enableGridGeometryCache>
54
61template<class GG>
63{
65 using GridView = typename GG::GridView;
66 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
67 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
68
69public:
71 using Element = typename GridView::template Codim<0>::Entity;
73 using SubControlVolume = typename GG::SubControlVolume;
75 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
77 using GridGeometry = GG;
78
80 static constexpr std::size_t maxNumElementScvs = 1;
82 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
83
86 : gridGeometryPtr_(&gridGeometry) {}
87
90 const SubControlVolume& scv(GridIndexType scvIdx) const
91 {
92 return gridGeometry().scv(scvIdx);
93 }
94
97 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
98 {
99 return gridGeometry().scvf(scvfIdx);
100 }
101
104 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
105 {
106 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
107 }
108
114 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
115 scvs(const CCTpfaFVElementGeometry& fvGeometry)
116 {
118 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
119 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
120 }
121
127 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
129 {
130 const auto& g = fvGeometry.gridGeometry();
131 const auto scvIdx = fvGeometry.scvIndices_[0];
133 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
134 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
135 }
136
138 std::size_t numScv() const
139 {
140 return scvIndices_.size();
141 }
142
144 std::size_t numScvf() const
145 {
146 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
147 }
148
155 {
156 this->bindElement(element);
157 return std::move(*this);
158 }
159
160 void bind(const Element& element) &
161 {
162 this->bindElement(element);
163 }
164
171 {
172 this->bindElement(element);
173 return std::move(*this);
174 }
175
177 void bindElement(const Element& element) &
178 {
179 element_ = element;
180 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
181 }
182
184 bool isBound() const
185 { return static_cast<bool>(element_); }
186
188 const Element& element() const
189 { return *element_; }
190
193 { return *gridGeometryPtr_; }
194
196 bool hasBoundaryScvf() const
197 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
198
199 typename Element::Geometry geometry(const SubControlVolume& scv) const
200 { return gridGeometryPtr_->element(scv.dofIndex()).geometry(); }
201
202 typename GridView::Intersection::Geometry geometry(const SubControlVolumeFace& scvf) const
203 {
204 const auto element = gridGeometryPtr_->element(scvf.insideScvIdx());
205 const auto& scvfIndices = gridGeometryPtr_->scvfIndicesOfScv(scvf.insideScvIdx());
206 const LocalIndexType localScvfIdx = Detail::Tpfa::findLocalIndex(scvf.index(), scvfIndices);
207 LocalIndexType localIdx = 0;
208 for (const auto& intersection : intersections(gridGeometryPtr_->gridView(), element))
209 {
210 if (intersection.neighbor() || intersection.boundary())
211 {
212 if (localIdx == localScvfIdx)
213 return intersection.geometry();
214 else
215 ++localIdx;
216 }
217 }
218
219 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf geometry");
220 }
221
222private:
223
224 std::optional<Element> element_;
225 std::array<GridIndexType, 1> scvIndices_;
226 const GridGeometry* gridGeometryPtr_;
227};
228
234template<class GG>
236{
238 using GridView = typename GG::GridView;
239 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
240 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
241
242 static const int dim = GridView::dimension;
243 static const int dimWorld = GridView::dimensionworld;
244
245public:
247 using Element = typename GridView::template Codim<0>::Entity;
249 using SubControlVolume = typename GG::SubControlVolume;
251 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
253 using GridGeometry = GG;
255 static constexpr std::size_t maxNumElementScvs = 1;
257 static constexpr std::size_t maxNumElementScvfs = 2*dim;
258
261 : gridGeometryPtr_(&gridGeometry) {}
262
265 const SubControlVolume& scv(GridIndexType scvIdx) const
266 {
267 if (scvIdx == scvIndices_[0])
268 return scvs_[0];
269 else
270 return neighborScvs_[Detail::Tpfa::findLocalIndex(scvIdx, neighborScvIndices_)];
271 }
272
275 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
276 {
277 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
278 if (it != scvfIndices_.end())
279 return scvfs_[std::distance(scvfIndices_.begin(), it)];
280 else
281 return neighborScvfs_[Detail::Tpfa::findLocalIndex(scvfIdx, neighborScvfIndices_)];
282 }
283
286 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
287 {
288 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
289 if (it != scvfIndices_.end())
290 {
291 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
292 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
293 }
294 else
295 {
296 const auto localScvfIdx = Detail::Tpfa::findLocalIndex(scvfIdx, neighborScvfIndices_);
297 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
298 if (localFlippedIndex < scvfs_.size())
299 return scvfs_[localFlippedIndex];
300 else
301 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
302 }
303 }
304
310 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
311 scvs(const ThisType& g)
312 {
313 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
314 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
315 }
316
322 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
323 scvfs(const ThisType& g)
324 {
325 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
326 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
327 }
328
330 std::size_t numScv() const
331 { return scvs_.size(); }
332
334 std::size_t numScvf() const
335 { return scvfs_.size(); }
336
343 {
344 this->bind_(element);
345 return std::move(*this);
346 }
347
348 void bind(const Element& element) &
349 {
350 this->bind_(element);
351 }
352
359 {
360 this->bindElement_(element);
361 return std::move(*this);
362 }
363
364 void bindElement(const Element& element) &
365 {
366 this->bindElement_(element);
367 }
368
370 bool isBound() const
371 { return static_cast<bool>(element_); }
372
374 const Element& element() const
375 { return *element_; }
376
379 { return *gridGeometryPtr_; }
380
382 bool hasBoundaryScvf() const
383 { return hasBoundaryScvf_; }
384
385 typename Element::Geometry geometry(const SubControlVolume& scv) const
386 { return gridGeometryPtr_->element(scv.dofIndex()).geometry(); }
387
388 typename GridView::Intersection::Geometry geometry(const SubControlVolumeFace& scvf) const
389 {
390 const auto element = gridGeometryPtr_->element(scvf.insideScvIdx());
391 const auto& scvfIndices = gridGeometryPtr_->scvfIndicesOfScv(scvf.insideScvIdx());
392 const LocalIndexType localScvfIdx = Detail::Tpfa::findLocalIndex(scvf.index(), scvfIndices);
393 LocalIndexType localIdx = 0;
394 for (const auto& intersection : intersections(gridGeometryPtr_->gridView(), element))
395 {
396 if (intersection.neighbor() || intersection.boundary())
397 {
398 if (localIdx == localScvfIdx)
399 return intersection.geometry();
400 else
401 ++localIdx;
402 }
403 }
404
405 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf geometry");
406 }
407
408private:
411 void bind_(const Element& element)
412 {
413 bindElement_(element);
414
415 neighborScvs_.reserve(element.subEntities(1));
416 neighborScvfIndices_.reserve(element.subEntities(1));
417 neighborScvfs_.reserve(element.subEntities(1));
418
419 std::vector<GridIndexType> handledNeighbors;
420 handledNeighbors.reserve(element.subEntities(1));
421
422 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
423 {
424 // for inner intersections and periodic (according to grid interface) intersections make neighbor geometry
425 if (intersection.neighbor())
426 {
427 const auto outside = intersection.outside();
428 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
429
430 // make outside geometries only if not done yet (could happen on non-conforming grids)
431 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
432 {
433 makeNeighborGeometries(outside, outsideIdx);
434 handledNeighbors.push_back(outsideIdx);
435 }
436 }
437 }
438
439 // build flip index set for network, surface, and periodic grids
440 if (dim < dimWorld || gridGeometry().isPeriodic())
441 {
442 flippedScvfIndices_.resize(scvfs_.size());
443 for (unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
444 {
445 const auto& scvf = scvfs_[localScvfIdx];
446 if (scvf.boundary())
447 continue;
448
449 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
450 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
451 {
452 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
453 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
454 {
455 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
456 {
457 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
458 break;
459 }
460 }
461 }
462 }
463
464 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
465 for (unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
466 {
467 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
468 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
469 for (unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
470 {
471 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
472 }
473 }
474 }
475 }
476
478 void bindElement_(const Element& element)
479 {
480 clear();
481 element_ = element;
482 scvfs_.reserve(element.subEntities(1));
483 scvfIndices_.reserve(element.subEntities(1));
484 makeElementGeometries(element);
485 }
486
487 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
488 {
489 for (unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
490 {
491 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
492 {
493 return scvfs_.size() + localNeighborScvfIdx;
494 }
495 }
496
497 // go over all potential scvfs of the outside scv
498 for (unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
499 {
500 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
501 for (unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
502 {
503 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
504 {
505 return localOutsideScvfIdx;
506 }
507 }
508 }
509
510 DUNE_THROW(Dune::InvalidStateException, "No flipped version of this scvf found!");
511 }
512
514 void makeElementGeometries(const Element& element)
515 {
516 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
517
518 const auto eIdx = gridGeometry().elementMapper().index(element);
519 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
520 scvIndices_[0] = eIdx;
521
522 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
523 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
524
525 // for network grids there might be multiple intersection with the same geometryInInside
526 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
527 // here we keep track of them
528 std::vector<bool> handledScvf;
529 if (dim < dimWorld)
530 handledScvf.resize(element.subEntities(1), false);
531
532 int scvfCounter = 0;
533 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
534 {
535 if (dim < dimWorld)
536 if (handledScvf[intersection.indexInInside()])
537 continue;
538
539 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
540 if (intersection.neighbor() || intersection.boundary())
541 {
542 ScvfGridIndexStorage scvIndices;
543 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
544 scvIndices[0] = eIdx;
545 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
546
547 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
548 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
549
550 scvfs_.emplace_back(intersection,
551 intersection.geometry(),
552 scvFaceIndices[scvfCounter],
553 scvIndices,
554 onBoundary);
555 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
556 scvfCounter++;
557
558 // for surface and network grids mark that we handled this face
559 if (dim < dimWorld)
560 handledScvf[intersection.indexInInside()] = true;
561 }
562 }
563 }
564
566 void makeNeighborGeometries(const Element& element, const GridIndexType eIdx)
567 {
568 using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
569
570 // create the neighbor scv
571 neighborScvs_.emplace_back(element.geometry(), eIdx);
572 neighborScvIndices_.push_back(eIdx);
573
574 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
575 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
576
577 // for network grids there might be multiple intersection with the same geometryInInside
578 // we identify those by the indexInInside for now (assumes conforming grids at branching facets)
579 // here we keep track of them
580 std::vector<bool> handledScvf;
581 if (dim < dimWorld)
582 handledScvf.resize(element.subEntities(1), false);
583
584 int scvfCounter = 0;
585 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
586 {
587 if (dim < dimWorld)
588 if (handledScvf[intersection.indexInInside()])
589 continue;
590
591 if (intersection.neighbor())
592 {
593 // this catches inner and periodic scvfs
594 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
595 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
596 {
597 // only create subcontrol faces where the outside element is the bound element
598 if (dim == dimWorld)
599 {
600 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
601 {
602 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
603 neighborScvfs_.emplace_back(intersection,
604 intersection.geometry(),
605 scvFaceIndices[scvfCounter],
606 scvIndices,
607 false);
608
609 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
610 }
611 }
612 // for network grids we can't use the intersection.outside() index as we can't assure that the
613 // first intersection with this indexInInside is the one that has our bound element as outside
614 // instead we check if the bound element's index is in the outsideScvIndices of the candidate scvf
615 // (will be optimized away for dim == dimWorld)
616 else
617 {
618 for (unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
619 {
620 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
621 {
622 ScvfGridIndexStorage scvIndices;
623 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
624 scvIndices[0] = eIdx;
625 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
626 neighborScvfs_.emplace_back(intersection,
627 intersection.geometry(),
628 scvFaceIndices[scvfCounter],
629 scvIndices,
630 false);
631
632 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
633 break;
634 }
635 }
636 }
637
638 // for surface and network grids mark that we handled this face
639 if (dim < dimWorld)
640 handledScvf[intersection.indexInInside()] = true;
641 scvfCounter++;
642 }
643 }
644
645 // only increase counter for boundary intersections
646 // (exclude periodic boundaries according to dune grid interface, they have been handled in neighbor==true)
647 else if (intersection.boundary())
648 {
649 if (dim < dimWorld)
650 handledScvf[intersection.indexInInside()] = true;
651 scvfCounter++;
652 }
653 }
654 }
655
657 void clear()
658 {
659 scvfIndices_.clear();
660 scvfs_.clear();
661 flippedScvfIndices_.clear();
662
663 neighborScvIndices_.clear();
664 neighborScvfIndices_.clear();
665 neighborScvs_.clear();
666 neighborScvfs_.clear();
667 flippedNeighborScvfIndices_.clear();
668
669 hasBoundaryScvf_ = false;
670 }
671
672 std::optional<Element> element_;
673 const GridGeometry* gridGeometryPtr_;
674
675 // local storage after binding an element
676 std::array<GridIndexType, 1> scvIndices_;
677 std::array<SubControlVolume, 1> scvs_;
678
679 std::vector<GridIndexType> scvfIndices_;
680 std::vector<SubControlVolumeFace> scvfs_;
681 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
682
683 std::vector<GridIndexType> neighborScvIndices_;
684 std::vector<SubControlVolume> neighborScvs_;
685
686 std::vector<GridIndexType> neighborScvfIndices_;
687 std::vector<SubControlVolumeFace> neighborScvfs_;
688 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
689
690 bool hasBoundaryScvf_ = false;
691};
692
693} // end namespace Dumux
694
695#endif
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:236
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:370
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:323
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:374
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:251
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:265
Element::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:385
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:330
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:260
GridView::Intersection::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:388
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:334
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:364
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:247
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:342
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:378
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:348
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:275
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:249
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:311
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:358
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:253
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:382
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:63
GridView::Intersection::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
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:170
Element::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:199
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:75
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:188
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:160
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:85
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:128
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:138
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:196
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:77
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:192
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:144
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:184
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:154
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:90
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:97
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:104
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:177
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:115
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:53
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:30
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:70
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:282
Defines the index types used for grid and local indices.
auto findLocalIndex(const GridIndexType idx, const std::vector< GridIndexType > &indices)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:33
Definition: adapt.hh:17
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28