version 3.8
discretization/pq1bubble/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_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <utility>
19
20#include <dune/common/exceptions.hh>
21#include <dune/geometry/type.hh>
22#include <dune/localfunctions/lagrange/pqkfactory.hh>
23
26
28
29namespace Dumux {
30
39template<class GG, bool enableGridGeometryCache>
41
43template<class GG>
45{
46 using GridView = typename GG::GridView;
47 static constexpr int dim = GridView::dimension;
48 static constexpr int dimWorld = GridView::dimensionworld;
49 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
50 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
51 using CoordScalar = typename GridView::ctype;
52 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
53 using GGCache = typename GG::Cache;
54 using GeometryHelper = typename GGCache::GeometryHelper;
55public:
57 using Element = typename GridView::template Codim<0>::Entity;
59 using SubControlVolume = typename GG::SubControlVolume;
61 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
63 using GridGeometry = GG;
65 // ToDo get this from GG
66 static constexpr std::size_t maxNumElementScvs = (1<<dim) + 1;
67
69 PQ1BubbleFVElementGeometry(const GGCache& ggCache)
70 : ggCache_(&ggCache)
71 {}
72
74 const SubControlVolume& scv(LocalIndexType scvIdx) const
75 {
76 return ggCache_->scvs(eIdx_)[scvIdx];
77 }
78
80 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
81 {
82 return ggCache_->scvfs(eIdx_)[scvfIdx];
83 }
84
90 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
92 {
93 using Iter = typename std::vector<SubControlVolume>::const_iterator;
94 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
95 return Dune::IteratorRange<Iter>(s.begin(), s.end());
96 }
97
103 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
105 {
106 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
107 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
108 return Dune::IteratorRange<Iter>(s.begin(), s.end());
109 }
110
112 const FeLocalBasis& feLocalBasis() const
113 {
114 return gridGeometry().feCache().get(element_->type()).localBasis();
115 }
116
118 std::size_t numScv() const
119 {
120 return ggCache_->scvs(eIdx_).size();
121 }
122
124 std::size_t numScvf() const
125 {
126 return ggCache_->scvfs(eIdx_).size();
127 }
128
135 {
136 this->bindElement(element);
137 return std::move(*this);
138 }
139
143 void bind(const Element& element) &
144 { this->bindElement(element); }
145
152 {
153 this->bindElement(element);
154 return std::move(*this);
155 }
156
160 void bindElement(const Element& element) &
161 {
162 element_ = element;
163 // cache element index
164 eIdx_ = 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 ggCache_->gridGeometry(); }
178
180 bool hasBoundaryScvf() const
181 { return ggCache_->hasBoundaryScvf(eIdx_); }
182
184 std::size_t elementIndex() const
185 { return eIdx_; }
186
188 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
189 {
190 if (scv.isOverlapping())
191 DUNE_THROW(Dune::NotImplemented, "Geometry of overlapping scv");
192
193 assert(isBound());
194 const auto geo = element().geometry();
195 const GeometryHelper helper(geo);
196 return {
197 helper.getScvGeometryType(scv.indexInElement()),
198 helper.getScvCorners(scv.indexInElement())
199 };
200 }
201
203 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
204 {
205 assert(isBound());
206 const auto geo = element().geometry();
207 if (scvf.boundary())
208 {
209 GeometryHelper helper(geo);
210 const auto localScvfIdx = scvf.index() - helper.numInteriorScvf();
211 const auto [localFacetIndex, isScvfLocalIdx]
212 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
213 return {
214 Dune::GeometryTypes::cube(dim-1),
215 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
216 };
217 }
218 else
219 {
220 GeometryHelper helper(geo);
221 return {
222 helper.getInteriorScvfGeometryType(scvf.index()),
223 helper.getScvfCorners(scvf.index())
224 };
225 }
226 }
227
228private:
229 const GGCache* ggCache_;
230 GridIndexType eIdx_;
231
232 std::optional<Element> element_;
233};
234
235} // end namespace Dumux
236
237#endif
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:91
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:104
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:57
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:180
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:74
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:63
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:112
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:61
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:59
PQ1BubbleFVElementGeometry 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/pq1bubble/fvelementgeometry.hh:151
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:160
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:143
PQ1BubbleFVElementGeometry 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/pq1bubble/fvelementgeometry.hh:134
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:203
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:118
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:188
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:80
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:168
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:184
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:124
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:176
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:172
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition: discretization/pq1bubble/fvelementgeometry.hh:40
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.
Defines the index types used for grid and local indices.
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