Loading [MathJax]/extensions/tex2jax.js
3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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;
73
76 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
77
80 : gridGeometryPtr_(&gridGeometry) {}
81
83 const SubControlVolume& scv(std::size_t scvIdx) const
84 { return gridGeometry().scvs(eIdx_)[scvIdx]; }
85
87 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
88 { return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
89
98 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
99 scvs(const BoxDfmFVElementGeometry& fvGeometry)
100 {
101 const auto& g = fvGeometry.gridGeometry();
102 using Iter = typename std::vector<SubControlVolume>::const_iterator;
103 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
104 }
105
114 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
116 {
117 const auto& g = fvGeometry.gridGeometry();
118 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
119 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
120 }
121
123 const FeLocalBasis& feLocalBasis() const
124 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
125
127 std::size_t numScv() const
128 { return gridGeometry().scvs(eIdx_).size(); }
129
131 std::size_t numScvf() const
132 { return gridGeometry().scvfs(eIdx_).size(); }
133
136 void bind(const Element& element)
137 {
138 this->bindElement(element);
139 }
140
147 void bindElement(const Element& element)
148 {
149 elemGeometryType_ = element.type();
150 eIdx_ = gridGeometry().elementMapper().index(element);
151 }
152
155 { return *gridGeometryPtr_; }
156
157private:
158 Dune::GeometryType elemGeometryType_;
159 const GridGeometry* gridGeometryPtr_;
160 GridIndexType eIdx_;
161};
162
164template<class GG>
166{
167 using GridView = typename GG::GridView;
168 static constexpr int dim = GridView::dimension;
169 static constexpr int dimWorld = GridView::dimensionworld;
170
171 using GridIndexType = typename GridView::IndexSet::IndexType;
172 using Element = typename GridView::template Codim<0>::Entity;
173
174 using CoordScalar = typename GridView::ctype;
175 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
176 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
177 using FaceReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim-1>;
178
179 using GeometryHelper = BoxDfmGeometryHelper<GridView, dim,
180 typename GG::SubControlVolume,
181 typename GG::SubControlVolumeFace>;
182public:
184 using SubControlVolume = typename GG::SubControlVolume;
186 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
188 using GridGeometry = GG;
191 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
192
195 : gridGeometryPtr_(&gridGeometry) {}
196
198 const SubControlVolume& scv(std::size_t scvIdx) const
199 { return scvs_[scvIdx]; }
200
202 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
203 { return scvfs_[scvfIdx]; }
204
213 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
214 scvs(const BoxDfmFVElementGeometry& fvGeometry)
215 {
216 using Iter = typename std::vector<SubControlVolume>::const_iterator;
217 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
218 }
219
228 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
230 {
231 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
232 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
233 }
234
236 const FeLocalBasis& feLocalBasis() const
237 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
238
240 std::size_t numScv() const
241 { return scvs_.size(); }
242
244 std::size_t numScvf() const
245 { return scvfs_.size(); }
246
249 void bind(const Element& element)
250 {
251 this->bindElement(element);
252 }
253
260 void bindElement(const Element& element)
261 {
262 eIdx_ = gridGeometry().elementMapper().index(element);
263 makeElementGeometries(element);
264 }
265
268 { return *gridGeometryPtr_; }
269
270private:
271
272 void makeElementGeometries(const Element& element)
273 {
274 auto eIdx = gridGeometry().elementMapper().index(element);
275
276 // get the element geometry
277 auto elementGeometry = element.geometry();
278 elemGeometryType_ = elementGeometry.type();
279 const auto referenceElement = ReferenceElements::general(elemGeometryType_);
280
281 // get the sub control volume geometries of this element
282 GeometryHelper geometryHelper(elementGeometry);
283
284 // construct the sub control volumes
285 scvs_.resize(elementGeometry.corners());
286 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
287 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
288 {
289 // get asssociated dof index
290 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
291
292 // add scv to the local container
293 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
294 scvLocalIdx,
295 eIdx,
296 dofIdxGlobal);
297 }
298
299 // construct the sub control volume faces
300 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
301 scvfs_.resize(numInnerScvf);
302
303 unsigned int scvfLocalIdx = 0;
304 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
305 {
306 // find the local scv indices this scvf is connected to
307 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
308 static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
309
310 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
311 element,
312 elementGeometry,
313 scvfLocalIdx,
314 std::move(localScvIndices));
315 }
316
317 // construct the ...
318 // ... sub-control volume faces on the domain boundary
319 // ... sub-control volumes on fracture facets
320 // ... sub-control volume faces on fracture facets
321 // NOTE We do not construct fracture scvfs on boundaries here!
322 // That means specifying Neumann fluxes on only the fractures is not possible
323 // However, it is difficult to interpret this here as a fracture ending on the boundary
324 // could also be connected to a facet which is both boundary and fracture at the same time!
325 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
326 // we would have to find only those fractures that are at the boundary and aren't connected
327 // to a fracture which is a boundary.
328 LocalIndexType scvLocalIdx = element.subEntities(dim);
329 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
330 {
331 // first, obtain all vertex indices on this intersection
332 const auto& isGeometry = intersection.geometry();
333 const auto numCorners = isGeometry.corners();
334 const auto idxInInside = intersection.indexInInside();
335
336 std::vector<GridIndexType> isVertexIndices(numCorners);
337 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
338 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
339 referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim),
340 dim);
341
342 if (intersection.boundary())
343 {
344 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
345 {
346 // find the scv this scvf is connected to
347 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
348 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
349
350 scvfs_.emplace_back(geometryHelper,
351 intersection,
352 isGeometry,
353 isScvfLocalIdx,
354 scvfLocalIdx,
355 std::move(localScvIndices));
356
357 // increment local counter
358 scvfLocalIdx++;
359 }
360 }
361
362 // maybe add fracture scvs & scvfs
363 if (this->gridGeometry().isOnFracture(element, intersection))
364 {
365 // add fracture scv for each vertex of intersection
366 const auto curNumScvs = scvs_.size();
367 scvs_.reserve(curNumScvs+numCorners);
368 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
369 scvs_.emplace_back(geometryHelper,
370 intersection,
371 isGeometry,
372 vIdxLocal,
373 static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
374 scvLocalIdx++,
375 idxInInside,
376 eIdx,
377 isVertexIndices[vIdxLocal]);
378
379 // add fracture scvf for each edge of the intersection in 3d
380 if (dim == 3)
381 {
382 const auto& faceReferenceElement = FaceReferenceElements::general(isGeometry.type());
383 for (unsigned int edgeIdx = 0; edgeIdx < faceReferenceElement.size(1); ++edgeIdx)
384 {
385 // inside/outside scv indices in face local node numbering
386 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceReferenceElement.subEntity(edgeIdx, 1, 0, dim-1)),
387 static_cast<LocalIndexType>(faceReferenceElement.subEntity(edgeIdx, 1, 1, dim-1))});
388
389 // add offset to get the right scv indices
390 std::for_each( localScvIndices.begin(),
391 localScvIndices.end(),
392 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
393
394 // add scvf
395 scvfs_.emplace_back(geometryHelper,
396 intersection,
397 isGeometry,
398 edgeIdx,
399 scvfLocalIdx++,
400 std::move(localScvIndices),
401 intersection.boundary());
402 }
403 }
404
405 // dim == 2, intersection is an edge, make 1 scvf
406 else
407 {
408 // inside/outside scv indices in face local node numbering
409 std::vector<LocalIndexType> localScvIndices({0, 1});
410
411 // add offset such that the fracture scvs above are addressed
412 std::for_each( localScvIndices.begin(),
413 localScvIndices.end(),
414 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
415
416 // add scvf
417 scvfs_.emplace_back(geometryHelper,
418 intersection,
419 isGeometry,
420 /*idxOnIntersection*/0,
421 scvfLocalIdx++,
422 std::move(localScvIndices),
423 intersection.boundary());
424 }
425 }
426 }
427 }
428
430 Dune::GeometryType elemGeometryType_;
431 GridIndexType eIdx_;
432
434 const GridGeometry* gridGeometryPtr_;
435
443 std::vector<SubControlVolume> scvs_;
444 std::vector<SubControlVolumeFace> scvfs_;
445};
446
447} // end namespace Dumux
448
449#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.
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:127
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:136
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:147
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:87
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:79
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:131
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:123
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:83
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:99
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:154
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:115
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:260
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:236
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:202
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:267
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:249
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:244
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:194
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:240
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:198
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:214
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:188
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:186
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:229
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:184
Create sub control volumes and sub control volume face geometries.
Definition: geometryhelper.hh:35