version 3.7
discretization/cellcentered/mpfa/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_CCMPFA_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <utility>
19
20#include <dune/common/exceptions.hh>
21#include <dune/common/iteratorrange.hh>
22
26
27namespace Dumux {
28
38template<class GG, bool enableGridGeometryCache>
40
47template<class GG>
49{
51 using GridView = typename GG::GridView;
52 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
53
54 static constexpr int dim = GridView::dimension;
55
56public:
58 using Element = typename GridView::template Codim<0>::Entity;
60 using SubControlVolume = typename GG::SubControlVolume;
62 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
64 using GridGeometry = GG;
66 static constexpr std::size_t maxNumElementScvs = 1;
68 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
69
72 : gridGeometryPtr_(&gridGeometry) {}
73
75 const SubControlVolume& scv(GridIndexType scvIdx) const
76 {
77 return gridGeometry().scv(scvIdx);
78 }
79
81 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
82 {
83 return gridGeometry().scvf(scvfIdx);
84 }
85
88 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
89 {
90 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
91 }
92
98 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
99 scvs(const CCMpfaFVElementGeometry& fvGeometry)
100 {
102 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
103 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
104 }
105
111 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
113 {
114 const auto& g = fvGeometry.gridGeometry();
115 const auto scvIdx = fvGeometry.scvIndices_[0];
117 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
118 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
119 }
120
122 std::size_t numScv() const
123 {
124 return scvIndices_.size();
125 }
126
128 std::size_t numScvf() const
129 {
130 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
131 }
132
139 {
140 this->bindElement(element);
141 return std::move(*this);
142 }
143
144 void bind(const Element& element) &
145 {
146 this->bindElement(element);
147 }
148
155 {
156 this->bindElement(element);
157 return std::move(*this);
158 }
159
161 void bindElement(const Element& element) &
162 {
163 element_ = element;
164 scvIndices_[0] = gridGeometry().elementMapper().index(element);
165 }
166
168 bool isBound() const
169 { return static_cast<bool>(element_); }
170
172 const Element& element() const
173 { return *element_; }
174
177 { return *gridGeometryPtr_; }
178
180 bool hasBoundaryScvf() const
181 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
182
183 // suppress warnings due to current implementation
184 // these interfaces should be used!
185 #pragma GCC diagnostic push
186 #pragma GCC diagnostic ignored "-Wdeprecated-declarations"
187
189 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
190 { return scv.geometry(); }
191
193 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
194 { return scvf.geometry(); }
195
196 #pragma GCC diagnostic pop
197
198private:
199
200 std::optional<Element> element_;
201 std::array<GridIndexType, 1> scvIndices_;
202 const GridGeometry* gridGeometryPtr_;
203};
204
210template<class GG>
212{
214 using GridView = typename GG::GridView;
215 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
216 using MpfaHelper = typename GG::MpfaHelper;
217
218 static const int dim = GridView::dimension;
219 static const int dimWorld = GridView::dimensionworld;
220 using CoordScalar = typename GridView::ctype;
221
222public:
224 using Element = typename GridView::template Codim<0>::Entity;
226 using SubControlVolume = typename GG::SubControlVolume;
228 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
230 using GridGeometry = GG;
232 static constexpr std::size_t maxNumElementScvs = 1;
234 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
235
238 : gridGeometryPtr_(&gridGeometry) {}
239
242 const SubControlVolume& scv(GridIndexType scvIdx) const
243 {
244 if (scvIdx == scvIndices_[0])
245 return scvs_[0];
246 else
247 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
248 }
249
252 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
253 {
254 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
255 if (it != scvfIndices_.end())
256 return scvfs_[std::distance(scvfIndices_.begin(), it)];
257 else
258 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
259 }
260
263 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
264 {
265 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
266 }
267
273 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
274 scvs(const ThisType& g)
275 {
276 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
277 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
278 }
279
285 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
286 scvfs(const ThisType& g)
287 {
288 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
289 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
290 }
291
293 std::size_t numScv() const
294 { return scvs_.size(); }
295
297 std::size_t numScvf() const
298 { return scvfs_.size(); }
299
306 {
307 this->bindElement_(element);
308 return std::move(*this);
309 }
310
311 void bindElement(const Element& element) &
312 {
313 this->bindElement_(element);
314 }
315
322 {
323 this->bind_(element);
324 return std::move(*this);
325 }
326
327 void bind(const Element& element) &
328 {
329 this->bind_(element);
330 }
331
333 bool isBound() const
334 { return static_cast<bool>(element_); }
335
337 const Element& element() const
338 { return *element_; }
339
342 { return *gridGeometryPtr_; }
343
345 bool hasBoundaryScvf() const
346 { return hasBoundaryScvf_; }
347
348 // suppress warnings due to current implementation
349 // these interfaces should be used!
350 #pragma GCC diagnostic push
351 #pragma GCC diagnostic ignored "-Wdeprecated-declarations"
352
354 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
355 { return scv.geometry(); }
356
358 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
359 { return scvf.geometry(); }
360
361 #pragma GCC diagnostic pop
362
363private:
364
367 void bind_(const Element& element)
368 {
369 // make inside geometries
370 bindElement_(element);
371
372 // get some references for convenience
373 const auto globalI = gridGeometry().elementMapper().index(element);
374 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
375
376 // reserve memory
377 const auto numNeighbors = assemblyMapI.size();
378 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
379 neighborScvs_.reserve(numNeighbors);
380 neighborScvIndices_.reserve(numNeighbors);
381 neighborScvfIndices_.reserve(numNeighborScvfs);
382 neighborScvfs_.reserve(numNeighborScvfs);
383
384 // make neighbor geometries
385 // use the assembly map to determine which faces are necessary
386 for (const auto& dataJ : assemblyMapI)
387 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
388 dataJ.globalJ,
389 dataJ.scvfsJ,
390 dataJ.additionalScvfs);
391
392 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
393 // //! on additional DOFs not included in the discretization schemes' occupation pattern
394 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
395 // if (!additionalDofDependencies.empty())
396 // {
397 // const auto newNumNeighbors = neighborScvs_.size() + additionalDofDependencies.size();
398 // neighborScvs_.reserve(newNumNeighbors);
399 // neighborScvIndices_.reserve(newNumNeighbors);
400 // for (auto globalJ : additionalDofDependencies)
401 // {
402 // neighborScvs_.emplace_back(gridGeometry().element(globalJ).geometry(), globalJ);
403 // neighborScvIndices_.emplace_back(globalJ);
404 // }
405 // }
406 }
407
409 void bindElement_(const Element& element)
410 {
411 clear();
412 element_ = element;
413 makeElementGeometries(element);
414 }
415
417 template<class DataJContainer>
418 std::size_t numNeighborScvfs_(const DataJContainer& dataJContainer)
419 {
420 std::size_t numNeighborScvfs = 0;
421 for (const auto& dataJ : dataJContainer)
422 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
423 return numNeighborScvfs;
424 }
425
427 void makeElementGeometries(const Element& element)
428 {
429 // make the scv
430 const auto eIdx = gridGeometry().elementMapper().index(element);
431 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
432 scvIndices_[0] = eIdx;
433
434 // get data on the scv faces
435 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
436 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
437
438 // the quadrature point parameterizaion to be used on scvfs
439 static const auto q = getParam<CoordScalar>("MPFA.Q");
440
441 // reserve memory for the scv faces
442 const auto numLocalScvf = scvFaceIndices.size();
443 scvfIndices_.reserve(numLocalScvf);
444 scvfs_.reserve(numLocalScvf);
445
446 // for network grids we only want to do one scvf per half facet
447 // this approach assumes conforming grids at branching facets
448 std::vector<bool> finishedFacets;
449 if (dim < dimWorld)
450 finishedFacets.resize(element.subEntities(1), false);
451
452 int scvfCounter = 0;
453 for (const auto& is : intersections(gridGeometry().gridView(), element))
454 {
455 // if we are dealing with a lower dimensional network
456 // only make a new scvf if we haven't handled it yet
457 if (dim < dimWorld)
458 {
459 const auto indexInInside = is.indexInInside();
460 if (finishedFacets[indexInInside])
461 continue;
462 else
463 finishedFacets[indexInInside] = true;
464 }
465
466 // if outside level > inside level, use the outside element in the following
467 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
468 const auto& e = useNeighbor ? is.outside() : element;
469 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
470 const auto eg = e.geometry();
471 const auto refElement = referenceElement(eg);
472
473 // Set up a container with all relevant positions for scvf corner computation
474 const auto numCorners = is.geometry().corners();
475 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
476 refElement,
477 indexInElement,
478 numCorners);
479
480 // make the scv faces belonging to each corner of the intersection
481 for (int c = 0; c < numCorners; ++c)
482 {
483 // get the global vertex index the scv face is connected to
484 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
485 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
486
487 // do not build scvfs connected to a processor boundary
488 if (gridGeometry().isGhostVertex(vIdxGlobal))
489 continue;
490
491 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
492
493 scvfs_.emplace_back(MpfaHelper(),
494 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
495 is,
496 vIdxGlobal,
497 vIdxLocal,
498 scvFaceIndices[scvfCounter],
499 eIdx,
500 neighborVolVarIndices[scvfCounter],
501 q,
502 is.boundary());
503
504 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
505 scvfCounter++;
506 }
507 }
508 }
509
511 template<typename IndexVector>
512 void makeNeighborGeometries(const Element& element,
513 GridIndexType eIdxGlobal,
514 const IndexVector& scvfIndices,
515 const IndexVector& additionalScvfs)
516 {
517 // create the neighbor scv if it doesn't exist yet
518 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
519 neighborScvIndices_.push_back(eIdxGlobal);
520
521 // get data on the scv faces
522 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
523 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
524
525 // the quadrature point parameterizaion to be used on scvfs
526 static const auto q = getParam<CoordScalar>("MPFA.Q");
527
528 // for network grids we only want to do one scvf per half facet
529 // this approach assumes conforming grids at branching facets
530 std::vector<bool> finishedFacets;
531 if (dim < dimWorld)
532 finishedFacets.resize(element.subEntities(1), false);
533
534 int scvfCounter = 0;
535 for (const auto& is : intersections(gridGeometry().gridView(), element))
536 {
537 // if we are dealing with a lower dimensional network
538 // only make a new scvf if we haven't handled it yet
539 if (dim < dimWorld)
540 {
541 auto indexInInside = is.indexInInside();
542 if(finishedFacets[indexInInside])
543 continue;
544 else
545 finishedFacets[indexInInside] = true;
546 }
547
548 // if outside level > inside level, use the outside element in the following
549 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
550 const auto& e = useNeighbor ? is.outside() : element;
551 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
552 const auto eg = e.geometry();
553 const auto refElement = referenceElement(eg);
554
555 // Set up a container with all relevant positions for scvf corner computation
556 const auto numCorners = is.geometry().corners();
557 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
558 refElement,
559 indexInElement,
560 numCorners);
561
562 // make the scv faces belonging to each corner of the intersection
563 for (int c = 0; c < numCorners; ++c)
564 {
565 // get the global vertex index the scv face is connected to
566 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
567 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
568
569 // do not build scvfs connected to a processor boundary
570 if (gridGeometry().isGhostVertex(vIdxGlobal))
571 continue;
572
573 // only build the scvf if it is in the list of necessary indices
574 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
575 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
576 {
577 // increment counter either way
578 scvfCounter++;
579 continue;
580 }
581
582 // build scvf
583 neighborScvfs_.emplace_back(MpfaHelper(),
584 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
585 is,
586 vIdxGlobal,
587 vIdxLocal,
588 scvFaceIndices[scvfCounter],
589 eIdxGlobal,
590 neighborVolVarIndices[scvfCounter],
591 q,
592 is.boundary());
593
594 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
595
596 // increment counter
597 scvfCounter++;
598 }
599 }
600 }
601
603 unsigned int findLocalIndex(const GridIndexType idx,
604 const std::vector<GridIndexType>& indices) const
605 {
606 auto it = std::find(indices.begin(), indices.end(), idx);
607 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
608 return std::distance(indices.begin(), it);
609 }
610
612 void clear()
613 {
614 scvfIndices_.clear();
615 scvfs_.clear();
616
617 neighborScvIndices_.clear();
618 neighborScvfIndices_.clear();
619 neighborScvs_.clear();
620 neighborScvfs_.clear();
621
622 hasBoundaryScvf_ = false;
623 }
624
625 const GridGeometry* gridGeometryPtr_;
626 std::optional<Element> element_;
627
628 // local storage after binding an element
629 std::array<GridIndexType, 1> scvIndices_;
630 std::vector<GridIndexType> scvfIndices_;
631 std::array<SubControlVolume, 1> scvs_;
632 std::vector<SubControlVolumeFace> scvfs_;
633
634 std::vector<GridIndexType> neighborScvIndices_;
635 std::vector<GridIndexType> neighborScvfIndices_;
636 std::vector<SubControlVolume> neighborScvs_;
637 std::vector<SubControlVolumeFace> neighborScvfs_;
638
639 bool hasBoundaryScvf_ = false;
640};
641
642} // end namespace
643
644#endif
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:212
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:252
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:242
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:286
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:354
void bindElement(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:311
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:341
CCMpfaFVElementGeometry 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/mpfa/fvelementgeometry.hh:321
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:263
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:226
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:224
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:293
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:228
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:337
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:333
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:327
CCMpfaFVElementGeometry 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/mpfa/fvelementgeometry.hh:305
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:297
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:345
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:358
GG GridGeometry
export type of finite volume grid geometries
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:230
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:274
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:237
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:49
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:75
CCMpfaFVElementGeometry 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/mpfa/fvelementgeometry.hh:154
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:193
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:161
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:128
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:71
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:172
CCMpfaFVElementGeometry 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/mpfa/fvelementgeometry.hh:138
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:168
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:99
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:144
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:180
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:62
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:112
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:81
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:88
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:176
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:122
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:58
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:60
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:189
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:64
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models This builds up th...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:39
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
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27