3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
28
29#include <optional>
30#include <utility>
31#include <unordered_map>
32#include <array>
33#include <vector>
34
35#include <dune/geometry/type.hh>
36#include <dune/localfunctions/lagrange/pqkfactory.hh>
37
41
42namespace Dumux {
43
52template<class GG, bool enableGridGeometryCache>
54
56template<class GG>
57class BoxFVElementGeometry<GG, true>
58{
59 using GridView = typename GG::GridView;
60 static constexpr int dim = GridView::dimension;
61 static constexpr int dimWorld = GridView::dimensionworld;
62 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
63 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
64 using CoordScalar = typename GridView::ctype;
65 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
66 using GGCache = typename GG::Cache;
67
68 using GeometryHelper = BoxGeometryHelper<GridView, dim,
69 typename GG::SubControlVolume,
70 typename GG::SubControlVolumeFace>;
71public:
73 using Element = typename GridView::template Codim<0>::Entity;
75 using SubControlVolume = typename GG::SubControlVolume;
77 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
79 using GridGeometry = GG;
81 static constexpr std::size_t maxNumElementScvs = (1<<dim);
82
87 BoxFVElementGeometry(const GGCache& ggCache)
88 : ggCache_(&ggCache) {}
89
90 [[deprecated("This Constructor is deprecated and will be removed after release 3.6. Always use localView(gridGeometry).")]]
92 : BoxFVElementGeometry(gg.cache_) {}
93
95 const SubControlVolume& scv(LocalIndexType scvIdx) const
96 {
97 return ggCache_->scvs(eIdx_)[scvIdx];
98 }
99
101 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
102 {
103 return ggCache_->scvfs(eIdx_)[scvfIdx];
104 }
105
111 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
112 scvs(const BoxFVElementGeometry& fvGeometry)
113 {
114 using Iter = typename std::vector<SubControlVolume>::const_iterator;
115 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
116 return Dune::IteratorRange<Iter>(s.begin(), s.end());
117 }
118
124 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
125 scvfs(const BoxFVElementGeometry& fvGeometry)
126 {
127 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
128 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
129 return Dune::IteratorRange<Iter>(s.begin(), s.end());
130 }
131
133 const FeLocalBasis& feLocalBasis() const
134 {
135 return gridGeometry().feCache().get(element_->type()).localBasis();
136 }
137
139 std::size_t numScv() const
140 {
141 return ggCache_->scvs(eIdx_).size();
142 }
143
145 std::size_t numScvf() const
146 {
147 return ggCache_->scvfs(eIdx_).size();
148 }
149
156 {
157 this->bindElement(element);
158 return std::move(*this);
159 }
160
164 void bind(const Element& element) &
165 { this->bindElement(element); }
166
173 {
174 this->bindElement(element);
175 return std::move(*this);
176 }
177
181 void bindElement(const Element& element) &
182 {
183 element_ = element;
184 // cache element index
185 eIdx_ = gridGeometry().elementMapper().index(element);
186 }
187
189 bool isBound() const
190 { return static_cast<bool>(element_); }
191
193 const Element& element() const
194 { return *element_; }
195
198 { return ggCache_->gridGeometry(); }
199
201 bool hasBoundaryScvf() const
202 { return ggCache_->hasBoundaryScvf(eIdx_); }
203
205 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
206 {
207 assert(isBound());
208 const auto geo = element().geometry();
209 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
210 }
211
213 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
214 {
215 assert(isBound());
216 const auto geo = element().geometry();
217 const GeometryHelper geometryHelper(geo);
218 if (scvf.boundary())
219 {
220 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
221 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
222 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
223 }
224 else
225 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
226 }
227
228private:
229 const GGCache* ggCache_;
230 GridIndexType eIdx_;
231
232 std::optional<Element> element_;
233};
234
236template<class GG>
237class BoxFVElementGeometry<GG, false>
238{
239 using GridView = typename GG::GridView;
240 static constexpr int dim = GridView::dimension;
241 static constexpr int dimWorld = GridView::dimensionworld;
242 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
243 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
244 using CoordScalar = typename GridView::ctype;
245 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
246 using GGCache = typename GG::Cache;
247
248 using GeometryHelper = BoxGeometryHelper<GridView, dim,
249 typename GG::SubControlVolume,
250 typename GG::SubControlVolumeFace>;
251public:
253 using Element = typename GridView::template Codim<0>::Entity;
255 using SubControlVolume = typename GG::SubControlVolume;
257 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
259 using GridGeometry = GG;
261 static constexpr std::size_t maxNumElementScvs = (1<<dim);
262
267 BoxFVElementGeometry(const GGCache& ggCache)
268 : ggCache_(&ggCache) {}
269
270 [[deprecated("This Constructor is deprecated and will be removed after release 3.6. Always use localView(gridGeometry).")]]
272 : BoxFVElementGeometry(gg.cache_) {}
273
275 const SubControlVolume& scv(LocalIndexType scvIdx) const
276 {
277 return scvs_[scvIdx];
278 }
279
281 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
282 {
283 return scvfs_[scvfIdx];
284 }
285
291 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
292 scvs(const BoxFVElementGeometry& fvGeometry)
293 {
294 using Iter = typename std::vector<SubControlVolume>::const_iterator;
295 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
296 }
297
303 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
304 scvfs(const BoxFVElementGeometry& fvGeometry)
305 {
306 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
307 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
308 }
309
311 const FeLocalBasis& feLocalBasis() const
312 {
313 return gridGeometry().feCache().get(element_->type()).localBasis();
314 }
315
317 std::size_t numScv() const
318 {
319 return scvs_.size();
320 }
321
323 std::size_t numScvf() const
324 {
325 return scvfs_.size();
326 }
327
334 {
335 this->bindElement(element);
336 return std::move(*this);
337 }
338
342 void bind(const Element& element) &
343 { this->bindElement(element); }
344
351 {
352 this->bindElement(element);
353 return std::move(*this);
354 }
355
359 void bindElement(const Element& element) &
360 {
361 element_ = element;
362 eIdx_ = gridGeometry().elementMapper().index(element);
363 makeElementGeometries_();
364 }
365
367 bool isBound() const
368 { return static_cast<bool>(element_); }
369
371 const Element& element() const
372 { return *element_; }
373
376 { return ggCache_->gridGeometry(); }
377
379 bool hasBoundaryScvf() const
380 { return hasBoundaryScvf_; }
381
383 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
384 {
385 assert(isBound());
386 const auto geo = element().geometry();
387 return { Dune::GeometryTypes::cube(dim), GeometryHelper(geo).getScvCorners(scv.indexInElement()) };
388 }
389
391 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
392 {
393 assert(isBound());
394 const auto geo = element().geometry();
395 const GeometryHelper geometryHelper(geo);
396 if (scvf.boundary())
397 {
398 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
399 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
400 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
401 }
402 else
403 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
404 }
405
406private:
407 void makeElementGeometries_()
408 {
409 hasBoundaryScvf_ = false;
410
411 // get the element geometry
412 const auto& element = *element_;
413 const auto elementGeometry = element.geometry();
414 const auto refElement = referenceElement(elementGeometry);
415
416 // get the sub control volume geometries of this element
417 GeometryHelper geometryHelper(elementGeometry);
418
419 // construct the sub control volumes
420 scvs_.resize(elementGeometry.corners());
421 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
422 {
423 // get associated dof index
424 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
425
426 // add scv to the local container
427 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
428 scvLocalIdx,
429 eIdx_,
430 dofIdxGlobal);
431 }
432
433 // construct the sub control volume faces
434 const auto numInnerScvf = geometryHelper.numInteriorScvf();
435 scvfs_.resize(numInnerScvf);
436 scvfBoundaryGeometryKeys_.clear();
437
438 LocalIndexType scvfLocalIdx = 0;
439 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
440 {
441 // find the local scv indices this scvf is connected to
442 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
443 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
444
445 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
446 element,
447 elementGeometry,
448 scvfLocalIdx,
449 std::move(localScvIndices),
450 false);
451 }
452
453 // construct the sub control volume faces on the domain boundary
454 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
455 {
456 if (intersection.boundary() && !intersection.neighbor())
457 {
458 const auto isGeometry = intersection.geometry();
459 hasBoundaryScvf_ = true;
460
461 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
462 {
463 // find the scv this scvf is connected to
464 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
465 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
466
467 scvfs_.emplace_back(geometryHelper,
468 intersection,
469 isGeometry,
470 isScvfLocalIdx,
471 scvfLocalIdx,
472 std::move(localScvIndices),
473 true);
474
475 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
476 static_cast<LocalIndexType>(intersection.indexInInside()),
477 static_cast<LocalIndexType>(isScvfLocalIdx)
478 }});
479
480 // increment local counter
481 scvfLocalIdx++;
482 }
483 }
484 }
485 }
486
488 GridIndexType eIdx_;
489 std::optional<Element> element_;
490
492 const GGCache* ggCache_;
493
495 std::vector<SubControlVolume> scvs_;
496 std::vector<SubControlVolumeFace> scvfs_;
497 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
498
499 bool hasBoundaryScvf_ = false;
500};
501
502} // end namespace Dumux
503
504#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.
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
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:53
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:189
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:133
BoxFVElementGeometry 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/box/fvelementgeometry.hh:172
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:193
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:205
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:112
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:139
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:95
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:145
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:75
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:197
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:201
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:125
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:164
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:213
BoxFVElementGeometry(const GridGeometry &gg)
Definition: discretization/box/fvelementgeometry.hh:91
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:77
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:181
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:101
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:87
BoxFVElementGeometry 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/box/fvelementgeometry.hh:155
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:73
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:79
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:367
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:323
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:259
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:391
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:379
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:292
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:317
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:359
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:281
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:267
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:383
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:375
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:253
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:304
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:311
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:371
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:275
BoxFVElementGeometry(const GridGeometry &gg)
Definition: discretization/box/fvelementgeometry.hh:271
BoxFVElementGeometry 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/box/fvelementgeometry.hh:333
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:257
BoxFVElementGeometry 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/box/fvelementgeometry.hh:350
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:255
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:342