3.5-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
224private:
225
226 const auto& scvfIndices_() const
227 {
228 return gridGeometry().scvfIndicesOfElement(eIdx_);
229 }
230
231 const Element* elementPtr_;
232 GridIndexType eIdx_;
233 const GridGeometry* gridGeometry_;
234};
235
241template<class GG>
242class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>
243{
244 using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>;
245
246 using GridView = typename GG::GridView;
247 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
248
249 //TODO include assert that checks for quad geometry
250
251 using StaticInfo = typename GG::StaticInformation;
252 static constexpr auto dim = StaticInfo::dim;
253 static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement;
254 static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv;
255 static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement;
256 static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement;
257 static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement;
258
259 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
260
261public:
263 using SubControlVolume = typename GG::SubControlVolume;
264 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
265 using Element = typename GridView::template Codim<0>::Entity;
266 using GridGeometry = GG;
267
269 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
270
272 : gridGeometry_(&gridGeometry)
273 , geometryHelper_(gridGeometry.gridView())
274 {}
275
277 const SubControlVolumeFace& scvf(const GridIndexType scvfIdx) const
278 { return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; }
279
285 friend inline auto
287 {
288 using Iter = typename decltype(fvGeometry.scvfs_)::const_iterator;
289 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
290 }
291
297 friend inline auto
299 {
300 using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>;
302 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
303 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
304 return Dune::IteratorRange<ScvfIterator>(begin, end);
305 }
306
308 const SubControlVolume& scv(const GridIndexType scvIdx) const
309 {
310 if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back())
311 return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)];
312 else
313 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
314 }
315
321 friend inline auto
323 {
324 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
325 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
326 }
327
330 {
331 assert(scvf.isLateral());
332 const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex());
333 return scvfs_[otherLocalIdx];
334 }
335
338 {
339 assert(scv.boundary());
340
341 // frontal boundary faces are always stored after the lateral faces
342 auto scvfIter = scvfs(*this, scv).begin();
343 const auto end = scvfs(*this, scv).end();
344 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
345 ++scvfIter;
346
347 assert(scvfIter->isFrontal());
348 assert(scvfIter->boundary());
349 return *scvfIter;
350 }
351
353 bool hasBoundaryScvf() const
354 { return gridGeometry().hasBoundaryScvf(eIdx_); }
355
357 std::size_t numScv() const
358 { return numScvsPerElement; }
359
361 std::size_t numScvf() const
362 { return scvfIndices_().size(); }
363
370 {
371 this->bind_(element);
372 return std::move(*this);
373 }
374
375 void bind(const Element& element) &
376 { this->bind_(element); }
377
384 {
385 typename GG::LocalIntersectionMapper localIsMapper;
386 localIsMapper.update(gridGeometry().gridView(), element);
387 this->bindElement_(element, localIsMapper);
388 return std::move(*this);
389 }
390
391 void bindElement(const Element& element) &
392 {
393 typename GG::LocalIntersectionMapper localIsMapper;
394 localIsMapper.update(gridGeometry().gridView(), element);
395 this->bindElement_(element, localIsMapper);
396 }
397
398 GridIndexType elementIndex() const
399 { return eIdx_; }
400
402 const Element& element() const
403 { return *elementPtr_; }
404
407 {
408 assert(gridGeometry_);
409 return *gridGeometry_;
410 }
411
414 { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); }
415
419 {
420 const SubControlVolumeFace& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
421 assert(!lateralOrthogonalScvf.boundary());
422
423 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
424
425 auto outsideFVGeometry = localView(gridGeometry());
426 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
427 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
428
429 for (const auto& otherScvf : scvfs(outsideFVGeometry))
430 if (otherScvf.localIndex() == otherLocalIdx)
431 return otherScvf;
432
433 DUNE_THROW(Dune::InvalidStateException, "No outside scvf found");
434 }
435
436private:
439 void bind_(const Element& element)
440 {
441 typename GG::LocalIntersectionMapper localIsMapper;
442 localIsMapper.update(gridGeometry().gridView(), element);
443
444 bindElement_(element, localIsMapper);
445 neighborScvIndices_.clear();
446 neighborScvs_.clear();
447
448 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
449 {
450 if (intersection.neighbor())
451 {
452 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
453 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
454 const auto& neighborElement = intersection.outside();
455 const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement);
456 const auto& neighborElementGeometry = neighborElement.geometry();
457
458 // todo: could be done easier?
459 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
460 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement);
461
462 typename GG::LocalIntersectionMapper localNeighborIsMapper;
463 localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement);
464
465 for (const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement))
466 {
467 const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside());
468 if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx)
469 {
470
471 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside());
472 neighborScvs_.push_back(SubControlVolume(
473 neighborElementGeometry,
474 neighborIntersection.geometry(),
475 globalScvIndicesOfNeighborElement[localNeighborScvIdx],
476 localNeighborScvIdx,
477 dofIndex,
478 Dumux::normalAxis(neighborIntersection.centerUnitOuterNormal()),
479 neighborElementIdx,
480 onDomainBoundary_(neighborIntersection)
481 ));
482
483 neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]);
484 }
485 }
486 }
487 }
488 }
489
491 void bindElement_(const Element& element, const typename GG::LocalIntersectionMapper& localIsMapper)
492 {
493 elementPtr_ = &element;
494 eIdx_ = gridGeometry().elementMapper().index(element);
495 std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement);
496 scvfs_.clear();
497 scvfs_.resize(minNumScvfsPerElement);
498
499 const auto& elementGeometry = element.geometry();
500
501 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
502 {
503 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
504 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
505 const auto& intersectionGeometry = intersection.geometry();
506 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
507
508 scvs_[localScvIdx] = SubControlVolume(
509 elementGeometry,
510 intersectionGeometry,
511 scvIndicesOfElement_[localScvIdx],
512 localScvIdx,
513 dofIndex,
514 Dumux::normalAxis(intersection.centerUnitOuterNormal()),
515 eIdx_,
516 onDomainBoundary_(intersection)
517 );
518
519 // the frontal sub control volume face at the element center
520 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
521 scvfs_[localScvfIdx] = SubControlVolumeFace(
522 elementGeometry,
523 intersectionGeometry,
524 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]},
525 localScvfIdx,
526 scvfIndices_()[localScvfIdx],
527 intersection.centerUnitOuterNormal(),
528 SubControlVolumeFace::FaceType::frontal,
529 SubControlVolumeFace::BoundaryType::interior
530 );
531 ++localScvfIdx;
532
533 // the lateral sub control volume faces
534 for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx),
535 [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
536 )
537 {
538 const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element);
539 const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element);
540 const auto& lateralFacetGeometry = lateralFacet.geometry();
541
542 // helper lambda to get the lateral scvf's global inside and outside scv indices
543 const auto globalScvIndicesForLateralFace = [&]
544 {
545 const auto globalOutsideScvIdx = [&]
546 {
547 if (lateralIntersection.neighbor())
548 {
549 const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside());
550 // todo: could be done easier?
551 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
552 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement);
553 return globalScvIndicesOfNeighborElement[localScvIdx];
554 }
555 else
556 return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]);
557 }();
558
559 return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx};
560 }();
561
562 const auto boundaryType = [&]
563 {
564 if (onProcessorBoundary_(lateralIntersection))
565 return SubControlVolumeFace::BoundaryType::processorBoundary;
566 else if (onDomainBoundary_(lateralIntersection))
567 return SubControlVolumeFace::BoundaryType::physicalBoundary;
568 else
569 return SubControlVolumeFace::BoundaryType::interior;
570 }();
571
572 scvfs_[localScvfIdx] = SubControlVolumeFace(
573 elementGeometry,
574 intersectionGeometry,
575 lateralFacetGeometry,
576 globalScvIndicesForLateralFace, // TODO higher order
577 localScvfIdx,
578 scvfIndices_()[localScvfIdx],
579 lateralIntersection.centerUnitOuterNormal(),
580 SubControlVolumeFace::FaceType::lateral,
581 boundaryType
582 );
583 ++localScvfIdx;
584 }
585 }
586
587 // do a second loop over all intersections to add frontal boundary faces
588 auto localScvfIdx = minNumScvfsPerElement;
589 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
590 {
591 // the frontal sub control volume face at a domain boundary (coincides with element face)
592 if (onDomainBoundary_(intersection))
593 {
594 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
595 // the frontal sub control volume face at the boundary
596 scvfs_.push_back(SubControlVolumeFace(
597 element.geometry(),
598 intersection.geometry(),
599 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]},
600 localScvfIdx,
601 scvfIndices_()[localScvfIdx],
602 intersection.centerUnitOuterNormal(),
603 SubControlVolumeFace::FaceType::frontal,
604 SubControlVolumeFace::BoundaryType::physicalBoundary
605 ));
606 ++localScvfIdx;
607 }
608 }
609
610 if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{})
611 {
612 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
613 if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement)
614 {
615 // make sure frontal boundary scvfs are sorted correctly
616 std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(),
617 [](const auto& scvfLeft, const auto& scvfRight)
618 { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); }
619 );
620 }
621 }
622 }
623
624 const auto& scvfIndices_() const
625 { return gridGeometry().scvfIndicesOfElement(eIdx_); }
626
627 template<class Entry, class Container>
628 const LocalIndexType findLocalIndex_(const Entry& entry,
629 const Container& container) const
630 {
631 auto it = std::find(container.begin(), container.end(), entry);
632 assert(it != container.end() && "Could not find the entry! Make sure to properly bind this class!");
633 return std::distance(container.begin(), it);
634 }
635
636 bool onDomainBoundary_(const typename GridView::Intersection& intersection) const
637 {
638 return !intersection.neighbor() && intersection.boundary();
639 }
640
641 bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const
642 {
643 return !intersection.neighbor() && !intersection.boundary();
644 }
645
646 Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_;
647
648 Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_;
649 Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_;
650
651 std::array<SubControlVolume, numScvsPerElement> scvs_;
652 std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_;
653
654 const GridGeometry* gridGeometry_;
655 const Element* elementPtr_;
656 GridIndexType eIdx_;
657 typename GridGeometry::GeometryHelper geometryHelper_;
658};
659
660} // end namespace Dumux
661
662#endif
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
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
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: adapt.hh:29
static std::size_t normalAxis(const Vector &v)
Definition: normalaxis.hh:39
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
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
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:243
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:361
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:353
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:286
void bindElement(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:391
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:263
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:369
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:413
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:383
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:266
GridIndexType elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:398
SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:418
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:265
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:329
const SubControlVolume & scv(const GridIndexType scvIdx) const
Get a half sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:308
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:406
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &g)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:322
void bind(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:375
const SubControlVolumeFace & scvf(const GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:277
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:402
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:337
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:298
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:271
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:357
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:264
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