3.1-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 <dune/common/version.hh>
33#include <dune/geometry/type.hh>
34#include <dune/geometry/referenceelements.hh>
35#include <dune/localfunctions/lagrange/pqkfactory.hh>
36
38#include "geometryhelper.hh"
39
40namespace Dumux {
41
51template<class GG, bool enableGridGeometryCache>
53
55template<class GG>
57{
58 using GridView = typename GG::GridView;
59 static constexpr int dim = GridView::dimension;
60 static constexpr int dimWorld = GridView::dimensionworld;
61 using GridIndexType = typename GridView::IndexSet::IndexType;
62 using Element = typename GridView::template Codim<0>::Entity;
63 using CoordScalar = typename GridView::ctype;
64 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
65 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
66public:
68 using SubControlVolume = typename GG::SubControlVolume;
70 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
72 using GridGeometry = GG;
74 using FVGridGeometry [[deprecated("Use GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry;
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(elemGeometryType_).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
138 void bind(const Element& element)
139 {
140 this->bindElement(element);
141 }
142
149 void bindElement(const Element& element)
150 {
151 elemGeometryType_ = element.type();
152 eIdx_ = gridGeometry().elementMapper().index(element);
153 }
154
156 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
158 { return *gridGeometryPtr_; }
159
162 { return *gridGeometryPtr_; }
163
164private:
165 Dune::GeometryType elemGeometryType_;
166 const GridGeometry* gridGeometryPtr_;
167 GridIndexType eIdx_;
168};
169
171template<class GG>
173{
174 using GridView = typename GG::GridView;
175 static constexpr int dim = GridView::dimension;
176 static constexpr int dimWorld = GridView::dimensionworld;
177
178 using GridIndexType = typename GridView::IndexSet::IndexType;
179 using Element = typename GridView::template Codim<0>::Entity;
180
181 using CoordScalar = typename GridView::ctype;
182 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
183 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
184 using FaceReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim-1>;
185
186 using GeometryHelper = BoxDfmGeometryHelper<GridView, dim,
187 typename GG::SubControlVolume,
188 typename GG::SubControlVolumeFace>;
189public:
191 using SubControlVolume = typename GG::SubControlVolume;
193 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
195 using GridGeometry = GG;
197 using FVGridGeometry [[deprecated("Use GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry;
200 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
201
204 : gridGeometryPtr_(&gridGeometry) {}
205
207 const SubControlVolume& scv(std::size_t scvIdx) const
208 { return scvs_[scvIdx]; }
209
211 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
212 { return scvfs_[scvfIdx]; }
213
222 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
223 scvs(const BoxDfmFVElementGeometry& fvGeometry)
224 {
225 using Iter = typename std::vector<SubControlVolume>::const_iterator;
226 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
227 }
228
237 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
239 {
240 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
241 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
242 }
243
245 const FeLocalBasis& feLocalBasis() const
246 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
247
249 std::size_t numScv() const
250 { return scvs_.size(); }
251
253 std::size_t numScvf() const
254 { return scvfs_.size(); }
255
258 void bind(const Element& element)
259 {
260 this->bindElement(element);
261 }
262
269 void bindElement(const Element& element)
270 {
271 eIdx_ = gridGeometry().elementMapper().index(element);
272 makeElementGeometries(element);
273 }
274
276 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
278 { return *gridGeometryPtr_; }
279
282 { return *gridGeometryPtr_; }
283
284private:
285
286 void makeElementGeometries(const Element& element)
287 {
288 auto eIdx = gridGeometry().elementMapper().index(element);
289
290 // get the element geometry
291 auto elementGeometry = element.geometry();
292 elemGeometryType_ = elementGeometry.type();
293 const auto referenceElement = ReferenceElements::general(elemGeometryType_);
294
295 // get the sub control volume geometries of this element
296 GeometryHelper geometryHelper(elementGeometry);
297
298 // construct the sub control volumes
299 scvs_.resize(elementGeometry.corners());
300 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
301 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
302 {
303 // get asssociated dof index
304 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
305
306 // add scv to the local container
307 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
308 scvLocalIdx,
309 eIdx,
310 dofIdxGlobal);
311 }
312
313 // construct the sub control volume faces
314 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
315 scvfs_.resize(numInnerScvf);
316
317 unsigned int scvfLocalIdx = 0;
318 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
319 {
320 // find the local scv indices this scvf is connected to
321 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
322 static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
323
324 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
325 element,
326 elementGeometry,
327 scvfLocalIdx,
328 std::move(localScvIndices));
329 }
330
331 // construct the ...
332 // ... sub-control volume faces on the domain boundary
333 // ... sub-control volumes on fracture facets
334 // ... sub-control volume faces on fracture facets
335 // NOTE We do not construct fracture scvfs on boundaries here!
336 // That means specifying Neumann fluxes on only the fractures is not possible
337 // However, it is difficult to interpret this here as a fracture ending on the boundary
338 // could also be connected to a facet which is both boundary and fracture at the same time!
339 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
340 // we would have to find only those fractures that are at the boundary and aren't connected
341 // to a fracture which is a boundary.
342 LocalIndexType scvLocalIdx = element.subEntities(dim);
343 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
344 {
345 // first, obtain all vertex indices on this intersection
346 const auto& isGeometry = intersection.geometry();
347 const auto numCorners = isGeometry.corners();
348 const auto idxInInside = intersection.indexInInside();
349
350 std::vector<GridIndexType> isVertexIndices(numCorners);
351 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
352 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
353 referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim),
354 dim);
355
356 if (intersection.boundary())
357 {
358 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
359 {
360 // find the scv this scvf is connected to
361 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
362 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
363
364 scvfs_.emplace_back(geometryHelper,
365 intersection,
366 isGeometry,
367 isScvfLocalIdx,
368 scvfLocalIdx,
369 std::move(localScvIndices));
370
371 // increment local counter
372 scvfLocalIdx++;
373 }
374 }
375
376 // maybe add fracture scvs & scvfs
377 if (this->gridGeometry().isOnFracture(element, intersection))
378 {
379 // add fracture scv for each vertex of intersection
380 const auto curNumScvs = scvs_.size();
381 scvs_.reserve(curNumScvs+numCorners);
382 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
383 scvs_.emplace_back(geometryHelper,
384 intersection,
385 isGeometry,
386 vIdxLocal,
387 static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
388 scvLocalIdx++,
389 idxInInside,
390 eIdx,
391 isVertexIndices[vIdxLocal]);
392
393 // add fracture scvf for each edge of the intersection in 3d
394 if (dim == 3)
395 {
396 const auto& faceReferenceElement = FaceReferenceElements::general(isGeometry.type());
397 for (unsigned int edgeIdx = 0; edgeIdx < faceReferenceElement.size(1); ++edgeIdx)
398 {
399 // inside/outside scv indices in face local node numbering
400 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceReferenceElement.subEntity(edgeIdx, 1, 0, dim-1)),
401 static_cast<LocalIndexType>(faceReferenceElement.subEntity(edgeIdx, 1, 1, dim-1))});
402
403 // add offset to get the right scv indices
404 std::for_each( localScvIndices.begin(),
405 localScvIndices.end(),
406 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
407
408 // add scvf
409 scvfs_.emplace_back(geometryHelper,
410 intersection,
411 isGeometry,
412 edgeIdx,
413 scvfLocalIdx++,
414 std::move(localScvIndices),
415 intersection.boundary());
416 }
417 }
418
419 // dim == 2, intersection is an edge, make 1 scvf
420 else
421 {
422 // inside/outside scv indices in face local node numbering
423 std::vector<LocalIndexType> localScvIndices({0, 1});
424
425 // add offset such that the fracture scvs above are addressed
426 std::for_each( localScvIndices.begin(),
427 localScvIndices.end(),
428 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
429
430 // add scvf
431 scvfs_.emplace_back(geometryHelper,
432 intersection,
433 isGeometry,
434 /*idxOnIntersection*/0,
435 scvfLocalIdx++,
436 std::move(localScvIndices),
437 intersection.boundary());
438 }
439 }
440 }
441 }
442
444 Dune::GeometryType elemGeometryType_;
445 GridIndexType eIdx_;
446
448 const GridGeometry* gridGeometryPtr_;
449
457 std::vector<SubControlVolume> scvs_;
458 std::vector<SubControlVolumeFace> scvfs_;
459};
460
461} // end namespace Dumux
462
463#endif
Class providing iterators over sub control volumes and sub control volume faces of an element.
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
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
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:52
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:70
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:129
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:138
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:149
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:72
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 GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:68
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:85
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:157
GridGeometry FVGridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:74
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:161
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
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:269
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:245
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:211
GridGeometry FVGridGeometry
export type of finite volume grid geometry
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:197
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:281
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:258
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:277
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:253
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:203
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:249
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:207
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:223
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:195
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:193
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:238
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:191
Create sub control volumes and sub control volume face geometries.
Definition: geometryhelper.hh:35