3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/box/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_FACETCOUPLING_BOX_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_FACETCOUPLING_BOX_FV_ELEMENT_GEOMETRY_HH
26
27#include <algorithm>
28
29#include <dune/geometry/type.hh>
30#include <dune/geometry/referenceelements.hh>
31
35
36namespace Dumux {
37
47template<class GG, bool enableGridGeometryCache>
49
51template<class GG>
53{
54 using GridView = typename GG::GridView;
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
57 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
58 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
59 using Element = typename GridView::template Codim<0>::Entity;
60 using CoordScalar = typename GridView::ctype;
61 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
62 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
63public:
65 using SubControlVolume = typename GG::SubControlVolume;
67 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
69 using GridGeometry = GG;
71 using FVGridGeometry [[deprecated("Use more general GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry;
73 static constexpr std::size_t maxNumElementScvs = (1<<dim);
74
77 : gridGeometryPtr_(&gridGeometry)
78 {}
79
81 const SubControlVolume& scv(LocalIndexType scvIdx) const
82 { return gridGeometry().scvs(eIdx_)[scvIdx]; }
83
85 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
86 { return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
87
93 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
95 {
96 const auto& g = fvGeometry.gridGeometry();
97 using Iter = typename std::vector<SubControlVolume>::const_iterator;
98 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
99 }
100
106 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
108 {
109 const auto& g = fvGeometry.gridGeometry();
110 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
111 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
112 }
113
115 const FeLocalBasis& feLocalBasis() const
116 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
117
119 std::size_t numScv() const
120 { return gridGeometry().scvs(eIdx_).size(); }
121
123 std::size_t numScvf() const
124 { return gridGeometry().scvfs(eIdx_).size(); }
125
129 void bind(const Element& element)
130 {
131 this->bindElement(element);
132 }
133
137 void bindElement(const Element& element)
138 {
139 elemGeometryType_ = element.type();
140 eIdx_ = gridGeometry().elementMapper().index(element);
141 }
142
144 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
146 { return *gridGeometryPtr_; }
147
150 { return *gridGeometryPtr_; }
151
152private:
153 Dune::GeometryType elemGeometryType_;
154 const GridGeometry* gridGeometryPtr_;
155
156 GridIndexType eIdx_;
157};
158
160template<class GG>
162{
163 using GridView = typename GG::GridView;
164 static constexpr int dim = GridView::dimension;
165 static constexpr int dimWorld = GridView::dimensionworld;
166
167 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
168 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
169 using Element = typename GridView::template Codim<0>::Entity;
170
171 using CoordScalar = typename GridView::ctype;
172 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
173 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
174
175 using GeometryHelper = BoxGeometryHelper<GridView, dim,
176 typename GG::SubControlVolume,
177 typename GG::SubControlVolumeFace>;
178public:
180 using SubControlVolume = typename GG::SubControlVolume;
182 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
184 using GridGeometry = GG;
186 using FVGridGeometry [[deprecated("Use more general GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry;
188 static constexpr std::size_t maxNumElementScvs = (1<<dim);
189
192 : gridGeometryPtr_(&gridGeometry)
193 {}
194
196 const SubControlVolume& scv(LocalIndexType scvIdx) const
197 { return scvs_[scvIdx]; }
198
200 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
201 { return scvfs_[scvfIdx]; }
202
208 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
210 {
211 using Iter = typename std::vector<SubControlVolume>::const_iterator;
212 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
213 }
214
220 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
222 {
223 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
224 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
225 }
226
228 const FeLocalBasis& feLocalBasis() const
229 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
230
232 std::size_t numScv() const
233 { return scvs_.size(); }
234
236 std::size_t numScvf() const
237 { return scvfs_.size(); }
238
242 void bind(const Element& element)
243 {
244 this->bindElement(element);
245 }
246
250 void bindElement(const Element& element)
251 {
252 eIdx_ = gridGeometry().elementMapper().index(element);
253 makeElementGeometries(element);
254 }
255
257 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
259 { return *gridGeometryPtr_; }
260
263 { return *gridGeometryPtr_; }
264
265private:
266
267 void makeElementGeometries(const Element& element)
268 {
269 auto eIdx = gridGeometry().elementMapper().index(element);
270
271 // get the element geometry
272 auto elementGeometry = element.geometry();
273 elemGeometryType_ = elementGeometry.type();
274 const auto referenceElement = ReferenceElements::general(elemGeometryType_);
275
276 // get the sub control volume geometries of this element
277 GeometryHelper geometryHelper(elementGeometry);
278
279 // construct the sub control volumes
280 scvs_.clear();
281 scvs_.reserve(elementGeometry.corners());
282 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
283 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
284 scvs_.emplace_back(geometryHelper,
285 scvLocalIdx,
286 eIdx,
287 gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim));
288
289 // construct the sub control volume faces
290 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
291 scvfs_.clear();
292 scvfs_.reserve(numInnerScvf);
293
294 unsigned int scvfLocalIdx = 0;
295 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
296 {
297 // find the local scv indices this scvf is connected to
298 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
299 static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
300
301 // create the sub-control volume face
302 scvfs_.emplace_back(geometryHelper,
303 element,
304 elementGeometry,
305 scvfLocalIdx,
306 std::move(localScvIndices));
307 }
308
309 // construct the sub control volume faces on the domain/interior boundaries
310 // skip handled facets (necessary for e.g. Dune::FoamGrid)
311 std::vector<unsigned int> handledFacets;
312 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
313 {
314 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
315 continue;
316
317 handledFacets.push_back(intersection.indexInInside());
318
319 // determine if all corners live on the facet grid
320 const auto isGeometry = intersection.geometry();
321 const auto numFaceCorners = isGeometry.corners();
322 const auto idxInInside = intersection.indexInInside();
323 const auto boundary = intersection.boundary();
324
325 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
326 for (int i = 0; i < numFaceCorners; ++i)
327 vIndicesLocal[i] = static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, i, dim));
328
329 // if all vertices are living on the facet grid, this is an interiour boundary
330 const bool isOnFacet = gridGeometry().isOnInteriorBoundary(element, intersection);
331
332 if (isOnFacet || boundary)
333 {
334 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
335 {
336 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
337 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
338
339 // create the sub-control volume face
340 scvfs_.emplace_back(geometryHelper,
341 intersection,
342 isGeometry,
343 isScvfLocalIdx,
344 scvfLocalIdx,
345 std::move(localScvIndices),
346 boundary,
347 isOnFacet);
348
349 // increment local counter
350 scvfLocalIdx++;
351 }
352 }
353 }
354 }
355
357 Dune::GeometryType elemGeometryType_;
358 GridIndexType eIdx_;
359
361 const GridGeometry* gridGeometryPtr_;
362
364 std::vector<SubControlVolume> scvs_;
365 std::vector<SubControlVolumeFace> scvfs_;
366};
367
368} // end namespace Dumux
369
370#endif
Defines the index types used for grid and local indices.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Class providing iterators over sub control volumes and sub control volume faces of an element.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:37
Base class for the element-local finite volume geometry for box models in the context of models consi...
Definition: multidomain/facet/box/fvelementgeometry.hh:48
GG GridGeometry
export type of finite volume grid geometry
Definition: multidomain/facet/box/fvelementgeometry.hh:69
void bind(const Element &element)
Definition: multidomain/facet/box/fvelementgeometry.hh:129
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvelementgeometry.hh:119
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: multidomain/facet/box/fvelementgeometry.hh:81
GridGeometry FVGridGeometry
export type of finite volume grid geometry
Definition: multidomain/facet/box/fvelementgeometry.hh:71
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFacetCouplingFVElementGeometry &fvGeometry)
Definition: multidomain/facet/box/fvelementgeometry.hh:107
void bindElement(const Element &element)
Definition: multidomain/facet/box/fvelementgeometry.hh:137
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: multidomain/facet/box/fvelementgeometry.hh:145
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: multidomain/facet/box/fvelementgeometry.hh:85
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: multidomain/facet/box/fvelementgeometry.hh:65
BoxFacetCouplingFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: multidomain/facet/box/fvelementgeometry.hh:76
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: multidomain/facet/box/fvelementgeometry.hh:149
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFacetCouplingFVElementGeometry &fvGeometry)
Definition: multidomain/facet/box/fvelementgeometry.hh:94
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: multidomain/facet/box/fvelementgeometry.hh:67
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: multidomain/facet/box/fvelementgeometry.hh:123
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: multidomain/facet/box/fvelementgeometry.hh:115
BoxFacetCouplingFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: multidomain/facet/box/fvelementgeometry.hh:191
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: multidomain/facet/box/fvelementgeometry.hh:196
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: multidomain/facet/box/fvelementgeometry.hh:182
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: multidomain/facet/box/fvelementgeometry.hh:228
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFacetCouplingFVElementGeometry &fvGeometry)
Definition: multidomain/facet/box/fvelementgeometry.hh:221
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvelementgeometry.hh:232
GridGeometry FVGridGeometry
export type of finite volume grid geometry
Definition: multidomain/facet/box/fvelementgeometry.hh:186
void bind(const Element &element)
Definition: multidomain/facet/box/fvelementgeometry.hh:242
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: multidomain/facet/box/fvelementgeometry.hh:236
void bindElement(const Element &element)
Definition: multidomain/facet/box/fvelementgeometry.hh:250
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: multidomain/facet/box/fvelementgeometry.hh:258
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFacetCouplingFVElementGeometry &fvGeometry)
Definition: multidomain/facet/box/fvelementgeometry.hh:209
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: multidomain/facet/box/fvelementgeometry.hh:180
GG GridGeometry
export type of finite volume grid geometry
Definition: multidomain/facet/box/fvelementgeometry.hh:184
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: multidomain/facet/box/fvelementgeometry.hh:262
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: multidomain/facet/box/fvelementgeometry.hh:200