3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/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 *****************************************************************************/
29#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH
30#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH
31
32#include <optional>
33#include <utility>
34
35#include <dune/common/version.hh>
36#include <dune/geometry/type.hh>
37#include <dune/localfunctions/lagrange/pqkfactory.hh>
38
40#include "geometryhelper.hh"
41
42namespace Dumux {
43
53template<class GG, bool enableGridGeometryCache>
55
57template<class GG>
59{
60 using GridView = typename GG::GridView;
61 static constexpr int dim = GridView::dimension;
62 static constexpr int dimWorld = GridView::dimensionworld;
63 using GridIndexType = typename GridView::IndexSet::IndexType;
64 using CoordScalar = typename GridView::ctype;
65 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
66public:
68 using Element = typename GridView::template Codim<0>::Entity;
70 using SubControlVolume = typename GG::SubControlVolume;
72 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
74 using GridGeometry = GG;
75
78 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
79
82 : gridGeometryPtr_(&gridGeometry) {}
83
85 const SubControlVolume& scv(std::size_t scvIdx) const
86 { return gridGeometry().scvs(eIdx_)[scvIdx]; }
87
89 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
90 { return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
91
100 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
101 scvs(const BoxDfmFVElementGeometry& fvGeometry)
102 {
103 const auto& g = fvGeometry.gridGeometry();
104 using Iter = typename std::vector<SubControlVolume>::const_iterator;
105 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
106 }
107
116 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
118 {
119 const auto& g = fvGeometry.gridGeometry();
120 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
121 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
122 }
123
125 const FeLocalBasis& feLocalBasis() const
126 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
127
129 std::size_t numScv() const
130 { return gridGeometry().scvs(eIdx_).size(); }
131
133 std::size_t numScvf() const
134 { return gridGeometry().scvfs(eIdx_).size(); }
135
142 {
143 this->bindElement(element);
144 return std::move(*this);
145 }
146
149 void bind(const Element& element) &
150 { this->bindElement(element); }
151
158 {
159 this->bindElement(element);
160 return std::move(*this);
161 }
162
169 void bindElement(const Element& element) &
170 {
171 element_ = element;
172 eIdx_ = gridGeometry().elementMapper().index(element);
173 }
174
177 { return *gridGeometryPtr_; }
178
180 bool isBound() const
181 { return static_cast<bool>(element_); }
182
184 const Element& element() const
185 { return *element_; }
186
187private:
188 const GridGeometry* gridGeometryPtr_;
189
190 std::optional<Element> element_;
191 GridIndexType eIdx_;
192};
193
195template<class GG>
197{
198 using GridView = typename GG::GridView;
199 static constexpr int dim = GridView::dimension;
200 static constexpr int dimWorld = GridView::dimensionworld;
201
202 using GridIndexType = typename GridView::IndexSet::IndexType;
203
204 using CoordScalar = typename GridView::ctype;
205 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
206 using GeometryHelper = BoxDfmGeometryHelper<GridView, dim,
207 typename GG::SubControlVolume,
208 typename GG::SubControlVolumeFace>;
209public:
211 using Element = typename GridView::template Codim<0>::Entity;
213 using SubControlVolume = typename GG::SubControlVolume;
215 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
217 using GridGeometry = GG;
220 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
221
224 : gridGeometryPtr_(&gridGeometry) {}
225
227 const SubControlVolume& scv(std::size_t scvIdx) const
228 { return scvs_[scvIdx]; }
229
231 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
232 { return scvfs_[scvfIdx]; }
233
242 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
243 scvs(const BoxDfmFVElementGeometry& fvGeometry)
244 {
245 using Iter = typename std::vector<SubControlVolume>::const_iterator;
246 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
247 }
248
257 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
259 {
260 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
261 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
262 }
263
265 const FeLocalBasis& feLocalBasis() const
266 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
267
269 std::size_t numScv() const
270 { return scvs_.size(); }
271
273 std::size_t numScvf() const
274 { return scvfs_.size(); }
275
282 {
283 this->bindElement(element);
284 return std::move(*this);
285 }
286
293 void bind(const Element& element) &
294 { this->bindElement(element); }
295
302 {
303 this->bindElement(element);
304 return std::move(*this);
305 }
306
311 void bindElement(const Element& element) &
312 {
313 element_ = element;
314 eIdx_ = gridGeometry().elementMapper().index(element);
315 makeElementGeometries_();
316 }
317
320 { return *gridGeometryPtr_; }
321
323 bool isBound() const
324 { return static_cast<bool>(element_); }
325
327 const Element& element() const
328 { return *element_; }
329
330private:
331
332 void makeElementGeometries_()
333 {
334 // get the element geometry
335 const auto& element = *element_;
336 const auto elementGeometry = element.geometry();
337 const auto refElement = referenceElement(element);
338
339 // get the sub control volume geometries of this element
340 GeometryHelper geometryHelper(elementGeometry);
341
342 // construct the sub control volumes
343 scvs_.resize(elementGeometry.corners());
344 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
345 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
346 {
347 // get asssociated dof index
348 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
349
350 // add scv to the local container
351 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
352 scvLocalIdx,
353 eIdx_,
354 dofIdxGlobal);
355 }
356
357 // construct the sub control volume faces
358 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
359 scvfs_.resize(numInnerScvf);
360
361 unsigned int scvfLocalIdx = 0;
362 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
363 {
364 // find the local scv indices this scvf is connected to
365 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
366 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
367
368 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
369 element,
370 elementGeometry,
371 scvfLocalIdx,
372 std::move(localScvIndices));
373 }
374
375 // construct the ...
376 // ... sub-control volume faces on the domain boundary
377 // ... sub-control volumes on fracture facets
378 // ... sub-control volume faces on fracture facets
379 // NOTE We do not construct fracture scvfs on boundaries here!
380 // That means specifying Neumann fluxes on only the fractures is not possible
381 // However, it is difficult to interpret this here as a fracture ending on the boundary
382 // could also be connected to a facet which is both boundary and fracture at the same time!
383 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
384 // we would have to find only those fractures that are at the boundary and aren't connected
385 // to a fracture which is a boundary.
386 LocalIndexType scvLocalIdx = element.subEntities(dim);
387 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
388 {
389 // first, obtain all vertex indices on this intersection
390 const auto& isGeometry = intersection.geometry();
391 const auto numCorners = isGeometry.corners();
392 const auto idxInInside = intersection.indexInInside();
393
394 std::vector<GridIndexType> isVertexIndices(numCorners);
395 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
396 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
397 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
398 dim);
399
400 if (intersection.boundary())
401 {
402 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
403 {
404 // find the scv this scvf is connected to
405 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
406 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
407
408 scvfs_.emplace_back(geometryHelper,
409 intersection,
410 isGeometry,
411 isScvfLocalIdx,
412 scvfLocalIdx,
413 std::move(localScvIndices));
414
415 // increment local counter
416 scvfLocalIdx++;
417 }
418 }
419
420 // maybe add fracture scvs & scvfs
421 if (this->gridGeometry().isOnFracture(element, intersection))
422 {
423 // add fracture scv for each vertex of intersection
424 const auto curNumScvs = scvs_.size();
425 scvs_.reserve(curNumScvs+numCorners);
426 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
427 scvs_.emplace_back(geometryHelper,
428 intersection,
429 isGeometry,
430 vIdxLocal,
431 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
432 scvLocalIdx++,
433 idxInInside,
434 eIdx_,
435 isVertexIndices[vIdxLocal]);
436
437 // add fracture scvf for each edge of the intersection in 3d
438 if (dim == 3)
439 {
440 const auto& faceRefElement = referenceElement(isGeometry);
441 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
442 {
443 // inside/outside scv indices in face local node numbering
444 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
445 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
446
447 // add offset to get the right scv indices
448 std::for_each( localScvIndices.begin(),
449 localScvIndices.end(),
450 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
451
452 // add scvf
453 scvfs_.emplace_back(geometryHelper,
454 intersection,
455 isGeometry,
456 edgeIdx,
457 scvfLocalIdx++,
458 std::move(localScvIndices),
459 intersection.boundary());
460 }
461 }
462
463 // dim == 2, intersection is an edge, make 1 scvf
464 else
465 {
466 // inside/outside scv indices in face local node numbering
467 std::vector<LocalIndexType> localScvIndices({0, 1});
468
469 // add offset such that the fracture scvs above are addressed
470 std::for_each( localScvIndices.begin(),
471 localScvIndices.end(),
472 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
473
474 // add scvf
475 scvfs_.emplace_back(geometryHelper,
476 intersection,
477 isGeometry,
478 /*idxOnIntersection*/0,
479 scvfLocalIdx++,
480 std::move(localScvIndices),
481 intersection.boundary());
482 }
483 }
484 }
485 }
486
488 std::optional<Element> element_;
489 GridIndexType eIdx_;
490
492 const GridGeometry* gridGeometryPtr_;
493
495 std::vector<SubControlVolume> scvs_;
496 std::vector<SubControlVolumeFace> scvfs_;
497};
498
499} // end namespace Dumux
500
501#endif
Class providing iterators over sub control volumes and sub control volume faces of an element.
Definition: adapt.hh:29
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:54
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:72
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:169
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:129
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:89
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:74
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:184
void bind(const Element &element) &
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:149
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:81
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:133
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:125
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:68
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:70
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:85
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:141
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:101
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:176
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:180
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:117
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:157
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:265
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:231
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:319
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:301
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:211
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:273
void bind(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:293
BoxDfmFVElementGeometry 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: porousmediumflow/boxdfm/fvelementgeometry.hh:281
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:327
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:223
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:269
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:323
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:227
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:243
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:217
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:215
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:258
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:213
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:311
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:35
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.