3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/fvgridgeometry.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_GRID_FVGEOMETRY_HH
30#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
31
32#include <unordered_map>
33
34#include <dune/localfunctions/lagrange/pqkfactory.hh>
35#include <dune/geometry/multilineargeometry.hh>
36#include <dune/grid/common/mcmgmapper.hh>
37
43
44#include "fvelementgeometry.hh"
45#include "geometryhelper.hh"
46#include "subcontrolvolume.hh"
48
49namespace Dumux {
50
59template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
61: public MapperTraits
62{
65
66 template<class GridGeometry, bool enableCache>
68
69 // Mapper type for mapping edges
70 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
71};
72
81template<class Scalar,
82 class GridView,
83 bool enableGridGeometryCache = false,
86
97template<class Scalar, class GV, class Traits>
98class BoxDfmFVGridGeometry<Scalar, GV, true, Traits>
99: public BaseGridGeometry<GV, Traits>
100{
103 using GridIndexType = typename GV::IndexSet::IndexType;
104
105 using Element = typename GV::template Codim<0>::Entity;
106 using CoordScalar = typename GV::ctype;
107 static const int dim = GV::dimension;
108 static const int dimWorld = GV::dimensionworld;
109 static_assert(dim == 2 || dim == 3, "The box-dfm GridGeometry is only implemented in 2 or 3 dimensions.");
110
111 using GeometryHelper = BoxDfmGeometryHelper<GV, dim,
112 typename Traits::SubControlVolume,
113 typename Traits::SubControlVolumeFace>;
114
115public:
117 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
118
120 using LocalView = typename Traits::template LocalView<ThisType, true>;
122 using SubControlVolume = typename Traits::SubControlVolume;
124 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
128 using DofMapper = typename Traits::VertexMapper;
130 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
132 using GridView = GV;
133
136 : ParentType(gridView) {}
137
140 const DofMapper& dofMapper() const
141 { return this->vertexMapper(); }
142
144 std::size_t numScv() const
145 { return numScv_; }
146
148 std::size_t numScvf() const
149 { return numScvf_; }
150
153 std::size_t numBoundaryScvf() const
154 { return numBoundaryScvf_; }
155
157 std::size_t numDofs() const
158 { return this->gridView().size(dim); }
159
161 template< class FractureGridAdapter >
162 void update(const FractureGridAdapter& fractureGridAdapter)
163 {
164 ParentType::update();
165
166 scvs_.clear();
167 scvfs_.clear();
168
169 auto numElements = this->gridView().size(0);
170 scvs_.resize(numElements);
171 scvfs_.resize(numElements);
172
173 boundaryDofIndices_.assign(numDofs(), false);
174 fractureDofIndices_.assign(this->gridView.size(dim), false);
175
176 numScv_ = 0;
177 numScvf_ = 0;
178 numBoundaryScvf_ = 0;
179 // Build the SCV and SCV faces
180 for (const auto& element : elements(this->gridView()))
181 {
182 // fill the element map with seeds
183 auto eIdx = this->elementMapper().index(element);
184
185 // count
186 numScv_ += element.subEntities(dim);
187 numScvf_ += element.subEntities(dim-1);
188
189 // get the element geometry
190 auto elementGeometry = element.geometry();
191 const auto refElement = referenceElement(elementGeometry);
192
193 // instantiate the geometry helper
194 GeometryHelper geometryHelper(elementGeometry);
195
196 // construct the sub control volumes
197 scvs_[eIdx].resize(elementGeometry.corners());
198 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
199 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
200 {
201 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
202
203 scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
204 scvLocalIdx,
205 eIdx,
206 dofIdxGlobal);
207 }
208
209 // construct the sub control volume faces
210 LocalIndexType scvfLocalIdx = 0;
211 scvfs_[eIdx].resize(element.subEntities(dim-1));
212 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
213 {
214 // find the global and local scv indices this scvf is belonging to
215 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
216 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
217
218 scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
219 element,
220 elementGeometry,
221 scvfLocalIdx,
222 std::move(localScvIndices));
223 }
224
225 // construct the ...
226 // ... sub-control volume faces on the domain boundary
227 // ... sub-control volumes on fracture facets
228 // ... sub-control volume faces on fracture facets
229 // NOTE We do not construct fracture scvfs on boundaries here!
230 // That means specifying Neumann fluxes on only the fractures is not possible
231 // However, it is difficult to interpret this here as a fracture ending on the boundary
232 // could also be connected to a facet which is both boundary and fracture at the same time!
233 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
234 // we would have to find only those fractures that are at the boundary and aren't connected
235 // to a fracture which is a boundary.
236 LocalIndexType scvLocalIdx = element.subEntities(dim);
237 for (const auto& intersection : intersections(this->gridView(), element))
238 {
239 // first, obtain all vertex indices on this intersection
240 const auto& isGeometry = intersection.geometry();
241 const auto numCorners = isGeometry.corners();
242 const auto idxInInside = intersection.indexInInside();
243
244 std::vector<GridIndexType> isVertexIndices(numCorners);
245 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
246 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
247 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
248 dim);
249 // maybe add boundary scvf
250 if (intersection.boundary() && !intersection.neighbor())
251 {
252 numScvf_ += isGeometry.corners();
253 numBoundaryScvf_ += isGeometry.corners();
254
255 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
256 {
257 // find the scvs this scvf is belonging to
258 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
259 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
260 scvfs_[eIdx].emplace_back(geometryHelper,
261 intersection,
262 isGeometry,
263 isScvfLocalIdx,
264 scvfLocalIdx++,
265 std::move(localScvIndices));
266 }
267
268 // add all vertices on the intersection to the set of
269 // boundary vertices
270 const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
271 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
272 {
273 const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
274 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
275 boundaryDofIndices_[vIdxGlobal] = true;
276 }
277 }
278 // make sure we have no periodic boundaries
279 else if (intersection.boundary() && intersection.neighbor())
280 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
281
282 // maybe add fracture scvs & scvfs
283 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
284 {
285 for (auto vIdx : isVertexIndices)
286 fractureDofIndices_[vIdx] = true;
287
288 // add fracture scv for each vertex of intersection
289 numScv_ += numCorners;
290 const auto curNumScvs = scvs_[eIdx].size();
291 scvs_[eIdx].reserve(curNumScvs+numCorners);
292 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
293 scvs_[eIdx].emplace_back(geometryHelper,
294 intersection,
295 isGeometry,
296 vIdxLocal,
297 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
298 scvLocalIdx++,
299 idxInInside,
300 eIdx,
301 isVertexIndices[vIdxLocal]);
302
303 // add fracture scvf for each edge of the intersection in 3d
304 if (dim == 3)
305 {
306 const auto& faceRefElement = referenceElement(isGeometry);
307 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
308 {
309 // inside/outside scv indices in face local node numbering
310 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
311 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
312
313 // add offset to get the right scv indices
314 std::for_each( localScvIndices.begin(),
315 localScvIndices.end(),
316 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
317
318 // add scvf
319 numScvf_++;
320 scvfs_[eIdx].emplace_back(geometryHelper,
321 intersection,
322 isGeometry,
323 edgeIdx,
324 scvfLocalIdx++,
325 std::move(localScvIndices),
326 intersection.boundary());
327 }
328 }
329
330 // dim == 2, intersection is an edge, make 1 scvf
331 else
332 {
333 // inside/outside scv indices in face local node numbering
334 std::vector<LocalIndexType> localScvIndices({0, 1});
335
336 // add offset such that the fracture scvs above are addressed
337 std::for_each( localScvIndices.begin(),
338 localScvIndices.end(),
339 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
340
341 // add scvf
342 numScvf_++;
343 scvfs_[eIdx].emplace_back(geometryHelper,
344 intersection,
345 isGeometry,
346 /*idxOnIntersection*/0,
347 scvfLocalIdx++,
348 std::move(localScvIndices),
349 intersection.boundary());
350 }
351 }
352 }
353 }
354 }
355
357 const FeCache& feCache() const { return feCache_; }
359 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const { return scvs_[eIdx]; }
361 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const { return scvfs_[eIdx]; }
363 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
365 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
367 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
368
370 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
371 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
372
374 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
375 { return std::unordered_map<std::size_t, std::size_t>(); }
376
377private:
378 const FeCache feCache_;
379
380 std::vector<std::vector<SubControlVolume>> scvs_;
381 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
382
383 // TODO do we need those?
384 std::size_t numScv_;
385 std::size_t numScvf_;
386 std::size_t numBoundaryScvf_;
387
388 // vertices on the boundary & fracture facets
389 std::vector<bool> boundaryDofIndices_;
390 std::vector<bool> fractureDofIndices_;
391};
392
400template<class Scalar, class GV, class Traits>
401class BoxDfmFVGridGeometry<Scalar, GV, false, Traits>
402: public BaseGridGeometry<GV, Traits>
403{
406 using GridIndexType = typename GV::IndexSet::IndexType;
407
408 static const int dim = GV::dimension;
409 static const int dimWorld = GV::dimensionworld;
410
411 using Element = typename GV::template Codim<0>::Entity;
412 using Intersection = typename GV::Intersection;
413 using CoordScalar = typename GV::ctype;
414
415public:
417 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
418
420 using LocalView = typename Traits::template LocalView<ThisType, false>;
422 using SubControlVolume = typename Traits::SubControlVolume;
424 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
428 using DofMapper = typename Traits::VertexMapper;
430 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
432 using GridView = GV;
433
436 : ParentType(gridView)
437 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
438 {}
439
442 const DofMapper& dofMapper() const
443 { return this->vertexMapper(); }
444
446 std::size_t numScv() const
447 { return numScv_; }
448
450 std::size_t numScvf() const
451 { return numScvf_; }
452
455 std::size_t numBoundaryScvf() const
456 { return numBoundaryScvf_; }
457
459 std::size_t numDofs() const
460 { return this->gridView().size(dim); }
461
463 template< class FractureGridAdapter >
464 void update(const FractureGridAdapter& fractureGridAdapter)
465 {
466 ParentType::update();
467
468 boundaryDofIndices_.assign(numDofs(), false);
469 fractureDofIndices_.assign(numDofs(), false);
470 facetOnFracture_.assign(this->gridView().size(1), false);
471
472 // save global data on the grid's scvs and scvfs
473 // TODO do we need those information?
474 numScv_ = 0;
475 numScvf_ = 0;
476 numBoundaryScvf_ = 0;
477 for (const auto& element : elements(this->gridView()))
478 {
479 numScv_ += element.subEntities(dim);
480 numScvf_ += element.subEntities(dim-1);
481
482 const auto elementGeometry = element.geometry();
483 const auto refElement = referenceElement(elementGeometry);
484
485 // store the sub control volume face indices on the domain boundary
486 for (const auto& intersection : intersections(this->gridView(), element))
487 {
488 // first, obtain all vertex indices on this intersection
489 const auto& isGeometry = intersection.geometry();
490 const auto numCorners = isGeometry.corners();
491 const auto idxInInside = intersection.indexInInside();
492
493 std::vector<GridIndexType> isVertexIndices(numCorners);
494 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
495 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
496 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
497 dim);
498
499 if (intersection.boundary() && !intersection.neighbor())
500 {
501 numScvf_ += numCorners;
502 numBoundaryScvf_ += numCorners;
503
504 // add all vertices on the intersection to the set of
505 // boundary vertices
506 const auto fIdx = intersection.indexInInside();
507 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
508 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
509 {
510 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
511 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
512 boundaryDofIndices_[vIdxGlobal] = true;
513 }
514 }
515
516 // make sure we have no periodic boundaries
517 else if (intersection.boundary() && intersection.neighbor())
518 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
519
520 // maybe add fracture scvs & scvfs
521 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
522 {
523 facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = true;
524 for (auto vIdx : isVertexIndices)
525 fractureDofIndices_[vIdx] = true;
526
527 const auto isGeometry = intersection.geometry();
528 numScv_ += isGeometry.corners();
529 numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
530 }
531 }
532 }
533 }
534
536 const FeCache& feCache() const { return feCache_; }
538 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
540 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
542 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
543
545 bool isOnFracture(const Element& element, const Intersection& intersection) const
546 { return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
547
549 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
550 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
551
553 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
554 { return std::unordered_map<std::size_t, std::size_t>(); }
555
556private:
557
558 const FeCache feCache_;
559
560 // Information on the global number of geometries
561 // TODO do we need those information?
562 std::size_t numScv_;
563 std::size_t numScvf_;
564 std::size_t numBoundaryScvf_;
565
566 // vertices on the boundary & fracture facets
567 std::vector<bool> boundaryDofIndices_;
568 std::vector<bool> fractureDofIndices_;
569
570 // facet mapper and markers which facets lie on interior boundaries
571 typename Traits::FacetMapper facetMapper_;
572 std::vector<bool> facetOnFracture_;
573};
574
575} // end namespace Dumux
576
577#endif
Defines the default element and vertex mapper types.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
The available discretization methods in Dumux.
Base class for grid geometries.
Helper classes to compute the integration elements.
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: common/pdesolver.hh:35
Base class for all finite volume grid geometries.
Definition: basegridgeometry.hh:49
GV GridView
export the grid view type
Definition: basegridgeometry.hh:64
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:51
The default traits for the box finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:62
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:70
Base class for the finite volume geometry vector for box schemes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:85
Base class for the finite volume geometry vector for box schemes that consider extra connectivity bet...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:100
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:122
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:374
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:361
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:157
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:128
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:135
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:130
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:359
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:126
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:120
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:140
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:148
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:144
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:370
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:153
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:363
void update(const FractureGridAdapter &fractureGridAdapter)
Update all fvElementGeometries (do this again after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:162
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:357
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:124
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:367
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:365
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:403
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:426
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:538
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:542
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:446
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:545
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:553
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:459
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:435
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:455
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:428
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:420
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:450
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:430
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:536
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:442
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:424
void update(const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (do this again after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:464
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:422
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:540
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:549
Create sub control volumes and sub control volume face geometries.
Definition: geometryhelper.hh:35
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:95
Class for a sub control volume face in the box discrete fracture method, i.e a part of the boundary o...
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:99
Base class for the local finite volume geometry for the box discrete fracture model.
the sub control volume for the box discrete fracture scheme
The sub control volume face class for the box discrete fracture model.