3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/facecentered/staggered/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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
26
27#include <utility>
28
29#include <dune/common/rangeutilities.hh>
30#include <dune/common/reservedvector.hh>
31
33#include <dune/common/iteratorrange.hh>
38
39namespace Dumux {
40
41template<class GG, bool cachingEnabled>
43
49template<class GG>
50class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true>
51{
52 using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true>;
53 using GridView = typename GG::GridView;
54 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
55 static constexpr auto numScvsPerElement = GG::StaticInformation::numScvsPerElement;
56
57public:
59 using SubControlVolume = typename GG::SubControlVolume;
60 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
61 using Element = typename GridView::template Codim<0>::Entity;
62 using GridGeometry = GG;
63
65 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
66
68 : gridGeometry_(&gridGeometry)
69 {}
70
72 const SubControlVolume& scv(GridIndexType scvIdx) const
73 { return gridGeometry().scv(scvIdx); }
74
76 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
77 { return gridGeometry().scvf(scvfIdx); }
78
81 {
82 assert(scvf.isLateral());
83 const auto otherGlobalIdx = scvfIndices_()[GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex())];
84 return gridGeometry().scvf(otherGlobalIdx);
85
86 }
87
90 {
91 assert(scv.boundary());
92
93 // frontal boundary faces are always stored after the lateral faces
94 auto scvfIter = scvfs(*this, scv).begin();
95 const auto end = scvfs(*this, scv).end();
96 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
97 ++scvfIter;
98
99 assert(scvfIter->isFrontal());
100 assert(scvfIter->boundary());
101 return *scvfIter;
102 }
103
109 friend inline auto
111 { return fvGeometry.gridGeometry().scvs(fvGeometry); }
112
118 friend inline auto
120 {
121 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
123 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(fvGeometry.scvfIndices_().begin(), fvGeometry),
124 ScvfIterator(fvGeometry.scvfIndices_().end(), fvGeometry));
125 }
126
132 friend inline auto
134 {
135 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
137 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
138 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
139 return Dune::IteratorRange<ScvfIterator>(begin, end);
140 }
141
143 std::size_t numScv() const
144 { return numScvsPerElement; }
145
147 std::size_t numScvf() const
148 { return scvfIndices_().size(); }
149
151 bool hasBoundaryScvf() const
152 { return gridGeometry().hasBoundaryScvf(eIdx_); }
153
160 {
161 this->bind(element);
162 return std::move(*this);
163 }
164
166 void bind(const Element& element) &
167 { this->bindElement(element); }
168
175 {
176 this->bindElement(element);
177 return std::move(*this);
178 }
179
181 void bindElement(const Element& element) &
182 {
183 elementPtr_ = &element;
184 eIdx_ = gridGeometry().elementMapper().index(element);
185 }
186
189 {
190 assert(gridGeometry_);
191 return *gridGeometry_;
192 }
193
194 std::size_t elementIndex() const
195 { return eIdx_; }
196
198 const Element& element() const
199 { return *elementPtr_; }
200
203 { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); }
204
207 {
208 const auto& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
209 assert(!lateralOrthogonalScvf.boundary());
210
211 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
212
213 auto outsideFVGeometry = localView(gridGeometry());
214 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
215 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
216
217 for (const auto& otherScvf : scvfs(outsideFVGeometry))
218 if (otherScvf.localIndex() == otherLocalIdx)
219 return otherScvf;
220
221 DUNE_THROW(Dune::InvalidStateException, "No outside scvf found");
222 }
223
224 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
225 { return scv.geometry(); }
226
227 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
228 { return scvf.geometry(); }
229
230private:
231
232 const auto& scvfIndices_() const
233 {
234 return gridGeometry().scvfIndicesOfElement(eIdx_);
235 }
236
237 const Element* elementPtr_;
238 GridIndexType eIdx_;
239 const GridGeometry* gridGeometry_;
240};
241
247template<class GG>
248class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>
249{
250 using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>;
251
252 using GridView = typename GG::GridView;
253 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
254
255 //TODO include assert that checks for quad geometry
256
257 using StaticInfo = typename GG::StaticInformation;
258 static constexpr auto dim = StaticInfo::dim;
259 static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement;
260 static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv;
261 static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement;
262 static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement;
263 static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement;
264
265 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
266
267public:
269 using SubControlVolume = typename GG::SubControlVolume;
270 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
271 using Element = typename GridView::template Codim<0>::Entity;
272 using GridGeometry = GG;
273
275 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
276
278 : gridGeometry_(&gridGeometry)
279 , geometryHelper_(gridGeometry.gridView())
280 {}
281
283 const SubControlVolumeFace& scvf(const GridIndexType scvfIdx) const
284 { return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; }
285
291 friend inline auto
293 {
294 using Iter = typename decltype(fvGeometry.scvfs_)::const_iterator;
295 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
296 }
297
303 friend inline auto
305 {
306 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
308 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
309 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
310 return Dune::IteratorRange<ScvfIterator>(begin, end);
311 }
312
314 const SubControlVolume& scv(const GridIndexType scvIdx) const
315 {
316 if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back())
317 return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)];
318 else
319 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
320 }
321
327 friend inline auto
329 {
330 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
331 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
332 }
333
336 {
337 assert(scvf.isLateral());
338 const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex());
339 return scvfs_[otherLocalIdx];
340 }
341
344 {
345 assert(scv.boundary());
346
347 // frontal boundary faces are always stored after the lateral faces
348 auto scvfIter = scvfs(*this, scv).begin();
349 const auto end = scvfs(*this, scv).end();
350 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
351 ++scvfIter;
352
353 assert(scvfIter->isFrontal());
354 assert(scvfIter->boundary());
355 return *scvfIter;
356 }
357
359 bool hasBoundaryScvf() const
360 { return gridGeometry().hasBoundaryScvf(eIdx_); }
361
363 std::size_t numScv() const
364 { return numScvsPerElement; }
365
367 std::size_t numScvf() const
368 { return scvfIndices_().size(); }
369
376 {
377 this->bind_(element);
378 return std::move(*this);
379 }
380
381 void bind(const Element& element) &
382 { this->bind_(element); }
383
390 {
391 typename GG::LocalIntersectionMapper localIsMapper;
392 localIsMapper.update(gridGeometry().gridView(), element);
393 this->bindElement_(element, localIsMapper);
394 return std::move(*this);
395 }
396
397 void bindElement(const Element& element) &
398 {
399 typename GG::LocalIntersectionMapper localIsMapper;
400 localIsMapper.update(gridGeometry().gridView(), element);
401 this->bindElement_(element, localIsMapper);
402 }
403
404 GridIndexType elementIndex() const
405 { return eIdx_; }
406
408 const Element& element() const
409 { return *elementPtr_; }
410
413 {
414 assert(gridGeometry_);
415 return *gridGeometry_;
416 }
417
420 { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); }
421
425 {
426 const SubControlVolumeFace& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
427 assert(!lateralOrthogonalScvf.boundary());
428
429 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
430
431 auto outsideFVGeometry = localView(gridGeometry());
432 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
433 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
434
435 for (const auto& otherScvf : scvfs(outsideFVGeometry))
436 if (otherScvf.localIndex() == otherLocalIdx)
437 return otherScvf;
438
439 DUNE_THROW(Dune::InvalidStateException, "No outside scvf found");
440 }
441
442private:
445 void bind_(const Element& element)
446 {
447 typename GG::LocalIntersectionMapper localIsMapper;
448 localIsMapper.update(gridGeometry().gridView(), element);
449
450 bindElement_(element, localIsMapper);
451 neighborScvIndices_.clear();
452 neighborScvs_.clear();
453
454 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
455 {
456 if (intersection.neighbor())
457 {
458 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
459 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
460 const auto& neighborElement = intersection.outside();
461 const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement);
462 const auto& neighborElementGeometry = neighborElement.geometry();
463
464 // todo: could be done easier?
465 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
466 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement);
467
468 typename GG::LocalIntersectionMapper localNeighborIsMapper;
469 localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement);
470
471 for (const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement))
472 {
473 const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside());
474 if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx)
475 {
476
477 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside());
478 neighborScvs_.push_back(SubControlVolume(
479 neighborElementGeometry,
480 neighborIntersection.geometry(),
481 globalScvIndicesOfNeighborElement[localNeighborScvIdx],
482 localNeighborScvIdx,
483 dofIndex,
484 Dumux::normalAxis(neighborIntersection.centerUnitOuterNormal()),
485 neighborElementIdx,
486 onDomainBoundary_(neighborIntersection)
487 ));
488
489 neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]);
490 }
491 }
492 }
493 }
494 }
495
497 void bindElement_(const Element& element, const typename GG::LocalIntersectionMapper& localIsMapper)
498 {
499 elementPtr_ = &element;
500 eIdx_ = gridGeometry().elementMapper().index(element);
501 std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement);
502 scvfs_.clear();
503 scvfs_.resize(minNumScvfsPerElement);
504
505 const auto& elementGeometry = element.geometry();
506
507 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
508 {
509 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
510 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
511 const auto& intersectionGeometry = intersection.geometry();
512 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
513
514 scvs_[localScvIdx] = SubControlVolume(
515 elementGeometry,
516 intersectionGeometry,
517 scvIndicesOfElement_[localScvIdx],
518 localScvIdx,
519 dofIndex,
520 Dumux::normalAxis(intersection.centerUnitOuterNormal()),
521 eIdx_,
522 onDomainBoundary_(intersection)
523 );
524
525 // the frontal sub control volume face at the element center
526 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
527 scvfs_[localScvfIdx] = SubControlVolumeFace(
528 elementGeometry,
529 intersectionGeometry,
530 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]},
531 localScvfIdx,
532 scvfIndices_()[localScvfIdx],
533 intersection.centerUnitOuterNormal(),
534 SubControlVolumeFace::FaceType::frontal,
535 SubControlVolumeFace::BoundaryType::interior
536 );
537 ++localScvfIdx;
538
539 // the lateral sub control volume faces
540 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx),
541 [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
542 )
543 {
544 const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element);
545 const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element);
546 const auto& lateralFacetGeometry = lateralFacet.geometry();
547
548 // helper lambda to get the lateral scvf's global inside and outside scv indices
549 const auto globalScvIndicesForLateralFace = [&]
550 {
551 const auto globalOutsideScvIdx = [&]
552 {
553 if (lateralIntersection.neighbor())
554 {
555 const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside());
556 // todo: could be done easier?
557 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
558 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement);
559 return globalScvIndicesOfNeighborElement[localScvIdx];
560 }
561 else
562 return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]);
563 }();
564
565 return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx};
566 }();
567
568 const auto boundaryType = [&]
569 {
570 if (onProcessorBoundary_(lateralIntersection))
571 return SubControlVolumeFace::BoundaryType::processorBoundary;
572 else if (onDomainBoundary_(lateralIntersection))
573 return SubControlVolumeFace::BoundaryType::physicalBoundary;
574 else
575 return SubControlVolumeFace::BoundaryType::interior;
576 }();
577
578 scvfs_[localScvfIdx] = SubControlVolumeFace(
579 elementGeometry,
580 intersectionGeometry,
581 lateralFacetGeometry,
582 globalScvIndicesForLateralFace, // TODO higher order
583 localScvfIdx,
584 scvfIndices_()[localScvfIdx],
585 lateralIntersection.centerUnitOuterNormal(),
586 SubControlVolumeFace::FaceType::lateral,
587 boundaryType
588 );
589 ++localScvfIdx;
590 }
591 }
592
593 // do a second loop over all intersections to add frontal boundary faces
594 auto localScvfIdx = minNumScvfsPerElement;
595 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
596 {
597 // the frontal sub control volume face at a domain boundary (coincides with element face)
598 if (onDomainBoundary_(intersection))
599 {
600 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
601 // the frontal sub control volume face at the boundary
602 scvfs_.push_back(SubControlVolumeFace(
603 element.geometry(),
604 intersection.geometry(),
605 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]},
606 localScvfIdx,
607 scvfIndices_()[localScvfIdx],
608 intersection.centerUnitOuterNormal(),
609 SubControlVolumeFace::FaceType::frontal,
610 SubControlVolumeFace::BoundaryType::physicalBoundary
611 ));
612 ++localScvfIdx;
613 }
614 }
615
616 if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{})
617 {
618 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
619 if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement)
620 {
621 // make sure frontal boundary scvfs are sorted correctly
622 std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(),
623 [](const auto& scvfLeft, const auto& scvfRight)
624 { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); }
625 );
626 }
627 }
628 }
629
630 const auto& scvfIndices_() const
631 { return gridGeometry().scvfIndicesOfElement(eIdx_); }
632
633 template<class Entry, class Container>
634 const LocalIndexType findLocalIndex_(const Entry& entry,
635 const Container& container) const
636 {
637 auto it = std::find(container.begin(), container.end(), entry);
638 assert(it != container.end() && "Could not find the entry! Make sure to properly bind this class!");
639 return std::distance(container.begin(), it);
640 }
641
642 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
643 {
644 return !intersection.neighbor() && intersection.boundary();
645 }
646
647 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
648 {
649 return !intersection.neighbor() && !intersection.boundary();
650 }
651
652 Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_;
653
654 Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_;
655 Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_;
656
657 std::array<SubControlVolume, numScvsPerElement> scvs_;
658 std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_;
659
660 const GridGeometry* gridGeometry_;
661 const Element* elementPtr_;
662 GridIndexType eIdx_;
663 typename GridGeometry::GeometryHelper geometryHelper_;
664};
665
666} // end namespace Dumux
667
668#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
static std::size_t normalAxis(const Vector &v)
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: normalaxis.hh:38
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
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:42
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:51
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:59
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:159
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:174
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:143
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:198
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:188
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:119
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:76
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:227
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:151
void bind(const Element &element) &
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:166
const SubControlVolumeFace & outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Returns the the scvf of neighbor element with the same integration point and unit outer normal.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:206
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:110
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:72
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:60
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:67
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:62
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:202
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:224
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:61
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:133
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:80
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:147
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:89
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:181
std::size_t elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:194
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:249
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:367
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:359
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:292
void bindElement(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:397
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:269
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:375
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:419
FaceCenteredStaggeredFVElementGeometry 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/facecentered/staggered/fvelementgeometry.hh:389
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:272
GridIndexType elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:404
SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:424
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:271
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:335
const SubControlVolume & scv(const GridIndexType scvIdx) const
Get a half sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:314
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:412
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &g)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:328
void bind(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:381
const SubControlVolumeFace & scvf(const GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:283
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:408
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:343
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:304
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:277
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:363
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:270
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82
Iterators over sub control volume faces of an fv geometry and a given sub control volume.
Definition: scvandscvfiterators.hh:122