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