3.3.0
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/localfunctions/lagrange/pqkfactory.hh>
35
37#include "geometryhelper.hh"
38
39namespace Dumux {
40
50template<class GG, bool enableGridGeometryCache>
52
54template<class GG>
56{
57 using GridView = typename GG::GridView;
58 static constexpr int dim = GridView::dimension;
59 static constexpr int dimWorld = GridView::dimensionworld;
60 using GridIndexType = typename GridView::IndexSet::IndexType;
61 using Element = typename GridView::template Codim<0>::Entity;
62 using CoordScalar = typename GridView::ctype;
63 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
64public:
66 using SubControlVolume = typename GG::SubControlVolume;
68 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
70 using GridGeometry = GG;
71
74 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
75
78 : gridGeometryPtr_(&gridGeometry) {}
79
81 const SubControlVolume& scv(std::size_t scvIdx) const
82 { return gridGeometry().scvs(eIdx_)[scvIdx]; }
83
85 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
86 { return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
87
96 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
97 scvs(const BoxDfmFVElementGeometry& fvGeometry)
98 {
99 const auto& g = fvGeometry.gridGeometry();
100 using Iter = typename std::vector<SubControlVolume>::const_iterator;
101 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
102 }
103
112 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
114 {
115 const auto& g = fvGeometry.gridGeometry();
116 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
117 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
118 }
119
121 const FeLocalBasis& feLocalBasis() const
122 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
123
125 std::size_t numScv() const
126 { return gridGeometry().scvs(eIdx_).size(); }
127
129 std::size_t numScvf() const
130 { return gridGeometry().scvfs(eIdx_).size(); }
131
134 void bind(const Element& element)
135 {
136 this->bindElement(element);
137 }
138
145 void bindElement(const Element& element)
146 {
147 elemGeometryType_ = element.type();
148 eIdx_ = gridGeometry().elementMapper().index(element);
149 }
150
153 { return *gridGeometryPtr_; }
154
155private:
156 Dune::GeometryType elemGeometryType_;
157 const GridGeometry* gridGeometryPtr_;
158 GridIndexType eIdx_;
159};
160
162template<class GG>
164{
165 using GridView = typename GG::GridView;
166 static constexpr int dim = GridView::dimension;
167 static constexpr int dimWorld = GridView::dimensionworld;
168
169 using GridIndexType = typename GridView::IndexSet::IndexType;
170 using Element = typename GridView::template Codim<0>::Entity;
171
172 using CoordScalar = typename GridView::ctype;
173 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
174 using GeometryHelper = BoxDfmGeometryHelper<GridView, dim,
175 typename GG::SubControlVolume,
176 typename GG::SubControlVolumeFace>;
177public:
179 using SubControlVolume = typename GG::SubControlVolume;
181 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
183 using GridGeometry = GG;
186 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
187
190 : gridGeometryPtr_(&gridGeometry) {}
191
193 const SubControlVolume& scv(std::size_t scvIdx) const
194 { return scvs_[scvIdx]; }
195
197 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
198 { return scvfs_[scvfIdx]; }
199
208 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
209 scvs(const BoxDfmFVElementGeometry& fvGeometry)
210 {
211 using Iter = typename std::vector<SubControlVolume>::const_iterator;
212 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
213 }
214
223 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
225 {
226 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
227 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
228 }
229
231 const FeLocalBasis& feLocalBasis() const
232 { return gridGeometry().feCache().get(elemGeometryType_).localBasis(); }
233
235 std::size_t numScv() const
236 { return scvs_.size(); }
237
239 std::size_t numScvf() const
240 { return scvfs_.size(); }
241
244 void bind(const Element& element)
245 {
246 this->bindElement(element);
247 }
248
255 void bindElement(const Element& element)
256 {
257 eIdx_ = gridGeometry().elementMapper().index(element);
258 makeElementGeometries(element);
259 }
260
263 { return *gridGeometryPtr_; }
264
265private:
266
267 void makeElementGeometries(const Element& element)
268 {
269 auto eIdx = gridGeometry().elementMapper().index(element);
270
271 // get the element geometry
272 auto elementGeometry = element.geometry();
273 elemGeometryType_ = elementGeometry.type();
274 const auto refElement = referenceElement(element);
275
276 // get the sub control volume geometries of this element
277 GeometryHelper geometryHelper(elementGeometry);
278
279 // construct the sub control volumes
280 scvs_.resize(elementGeometry.corners());
281 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
282 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
283 {
284 // get asssociated dof index
285 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
286
287 // add scv to the local container
288 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
289 scvLocalIdx,
290 eIdx,
291 dofIdxGlobal);
292 }
293
294 // construct the sub control volume faces
295 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
296 scvfs_.resize(numInnerScvf);
297
298 unsigned int scvfLocalIdx = 0;
299 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
300 {
301 // find the local scv indices this scvf is connected to
302 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
303 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
304
305 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
306 element,
307 elementGeometry,
308 scvfLocalIdx,
309 std::move(localScvIndices));
310 }
311
312 // construct the ...
313 // ... sub-control volume faces on the domain boundary
314 // ... sub-control volumes on fracture facets
315 // ... sub-control volume faces on fracture facets
316 // NOTE We do not construct fracture scvfs on boundaries here!
317 // That means specifying Neumann fluxes on only the fractures is not possible
318 // However, it is difficult to interpret this here as a fracture ending on the boundary
319 // could also be connected to a facet which is both boundary and fracture at the same time!
320 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
321 // we would have to find only those fractures that are at the boundary and aren't connected
322 // to a fracture which is a boundary.
323 LocalIndexType scvLocalIdx = element.subEntities(dim);
324 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
325 {
326 // first, obtain all vertex indices on this intersection
327 const auto& isGeometry = intersection.geometry();
328 const auto numCorners = isGeometry.corners();
329 const auto idxInInside = intersection.indexInInside();
330
331 std::vector<GridIndexType> isVertexIndices(numCorners);
332 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
333 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
334 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
335 dim);
336
337 if (intersection.boundary())
338 {
339 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
340 {
341 // find the scv this scvf is connected to
342 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
343 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
344
345 scvfs_.emplace_back(geometryHelper,
346 intersection,
347 isGeometry,
348 isScvfLocalIdx,
349 scvfLocalIdx,
350 std::move(localScvIndices));
351
352 // increment local counter
353 scvfLocalIdx++;
354 }
355 }
356
357 // maybe add fracture scvs & scvfs
358 if (this->gridGeometry().isOnFracture(element, intersection))
359 {
360 // add fracture scv for each vertex of intersection
361 const auto curNumScvs = scvs_.size();
362 scvs_.reserve(curNumScvs+numCorners);
363 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
364 scvs_.emplace_back(geometryHelper,
365 intersection,
366 isGeometry,
367 vIdxLocal,
368 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
369 scvLocalIdx++,
370 idxInInside,
371 eIdx,
372 isVertexIndices[vIdxLocal]);
373
374 // add fracture scvf for each edge of the intersection in 3d
375 if (dim == 3)
376 {
377 const auto& faceRefElement = referenceElement(isGeometry);
378 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
379 {
380 // inside/outside scv indices in face local node numbering
381 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
382 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
383
384 // add offset to get the right scv indices
385 std::for_each( localScvIndices.begin(),
386 localScvIndices.end(),
387 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
388
389 // add scvf
390 scvfs_.emplace_back(geometryHelper,
391 intersection,
392 isGeometry,
393 edgeIdx,
394 scvfLocalIdx++,
395 std::move(localScvIndices),
396 intersection.boundary());
397 }
398 }
399
400 // dim == 2, intersection is an edge, make 1 scvf
401 else
402 {
403 // inside/outside scv indices in face local node numbering
404 std::vector<LocalIndexType> localScvIndices({0, 1});
405
406 // add offset such that the fracture scvs above are addressed
407 std::for_each( localScvIndices.begin(),
408 localScvIndices.end(),
409 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
410
411 // add scvf
412 scvfs_.emplace_back(geometryHelper,
413 intersection,
414 isGeometry,
415 /*idxOnIntersection*/0,
416 scvfLocalIdx++,
417 std::move(localScvIndices),
418 intersection.boundary());
419 }
420 }
421 }
422 }
423
425 Dune::GeometryType elemGeometryType_;
426 GridIndexType eIdx_;
427
429 const GridGeometry* gridGeometryPtr_;
430
438 std::vector<SubControlVolume> scvs_;
439 std::vector<SubControlVolumeFace> scvfs_;
440};
441
442} // end namespace Dumux
443
444#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:51
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:68
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:125
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:134
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:145
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:85
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:70
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:77
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:129
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:121
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:66
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:81
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:97
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:152
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:113
void bindElement(const Element &element)
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:255
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:231
const SubControlVolumeFace & scvf(std::size_t scvfIdx) const
Get a sub control volume face with a local scvf index.
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:262
void bind(const Element &element)
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:244
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:239
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:189
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:235
const SubControlVolume & scv(std::size_t scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:193
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:209
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:183
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:181
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:224
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:179
Create sub control volumes and sub control volume face geometries.
Definition: geometryhelper.hh:35