version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
17#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
18#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
19
20#include <utility>
21#include <unordered_map>
22
23#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
24#include <dune/geometry/multilineargeometry.hh>
25#include <dune/grid/common/mcmgmapper.hh>
26
32
33#include "fvelementgeometry.hh"
34#include "geometryhelper.hh"
35#include "subcontrolvolume.hh"
37
38namespace Dumux {
39
40namespace Detail {
41template<class GV, class T>
42using BoxDfmGeometryHelper_t = Dune::Std::detected_or_t<
45 T
46>;
47} // end namespace Detail
48
57template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
59: public MapperTraits
60{
63
64 template<class GridGeometry, bool enableCache>
66
67 // Mapper type for mapping edges
68 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
69};
70
79template<class Scalar,
80 class GridView,
81 bool enableGridGeometryCache = false,
84
95template<class Scalar, class GV, class Traits>
96class BoxDfmFVGridGeometry<Scalar, GV, true, Traits>
97: public BaseGridGeometry<GV, Traits>
98{
101 using GridIndexType = typename GV::IndexSet::IndexType;
102
103 using Element = typename GV::template Codim<0>::Entity;
104 using CoordScalar = typename GV::ctype;
105 static const int dim = GV::dimension;
106 static const int dimWorld = GV::dimensionworld;
107 static_assert(dim == 2 || dim == 3, "The box-dfm GridGeometry is only implemented in 2 or 3 dimensions.");
108
109public:
112 static constexpr DiscretizationMethod discMethod{};
113
115 using LocalView = typename Traits::template LocalView<ThisType, true>;
117 using SubControlVolume = typename Traits::SubControlVolume;
119 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
123 using DofMapper = typename Traits::VertexMapper;
125 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
127 using GridView = GV;
130
132 template< class FractureGridAdapter >
133 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
134 : ParentType(gridView)
135 {
136 update_(fractureGridAdapter);
137 }
138
141 const DofMapper& dofMapper() const
142 { return this->vertexMapper(); }
143
145 std::size_t numScv() const
146 { return numScv_; }
147
149 std::size_t numScvf() const
150 { return numScvf_; }
151
154 std::size_t numBoundaryScvf() const
155 { return numBoundaryScvf_; }
156
158 std::size_t numDofs() const
159 { return this->gridView().size(dim); }
160
162 template< class FractureGridAdapter >
163 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
164 {
165 ParentType::update(gridView);
166 update_(fractureGridAdapter);
167 }
168
170 template< class FractureGridAdapter >
171 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
172 {
173 ParentType::update(std::move(gridView));
174 update_(fractureGridAdapter);
175 }
176
178 const FeCache& feCache() const { return feCache_; }
180 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const { return scvs_[eIdx]; }
182 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const { return scvfs_[eIdx]; }
184 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
186 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
188 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
189
191 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
192 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
193
194private:
195
196 template< class FractureGridAdapter >
197 void update_(const FractureGridAdapter& fractureGridAdapter)
198 {
199 scvs_.clear();
200 scvfs_.clear();
201
202 auto numElements = this->gridView().size(0);
203 scvs_.resize(numElements);
204 scvfs_.resize(numElements);
205
206 boundaryDofIndices_.assign(numDofs(), false);
207 fractureDofIndices_.assign(this->gridView.size(dim), false);
208
209 numScv_ = 0;
210 numScvf_ = 0;
211 numBoundaryScvf_ = 0;
212 // Build the SCV and SCV faces
213 for (const auto& element : elements(this->gridView()))
214 {
215 // fill the element map with seeds
216 auto eIdx = this->elementMapper().index(element);
217
218 // count
219 numScv_ += element.subEntities(dim);
220 numScvf_ += element.subEntities(dim-1);
221
222 // get the element geometry
223 auto elementGeometry = element.geometry();
224 const auto refElement = referenceElement(elementGeometry);
225
226 // instantiate the geometry helper
227 GeometryHelper geometryHelper(elementGeometry);
228
229 // construct the sub control volumes
230 scvs_[eIdx].resize(elementGeometry.corners());
231 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
232 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
233 {
234 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
235
236 scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
237 scvLocalIdx,
238 eIdx,
239 dofIdxGlobal);
240 }
241
242 // construct the sub control volume faces
243 LocalIndexType scvfLocalIdx = 0;
244 scvfs_[eIdx].resize(element.subEntities(dim-1));
245 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
246 {
247 // find the global and local scv indices this scvf is belonging to
248 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
249 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
250
251 scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
252 element,
253 elementGeometry,
254 scvfLocalIdx,
255 std::move(localScvIndices));
256 }
257
258 // construct the ...
259 // ... sub-control volume faces on the domain boundary
260 // ... sub-control volumes on fracture facets
261 // ... sub-control volume faces on fracture facets
262 // NOTE We do not construct fracture scvfs on boundaries here!
263 // That means specifying Neumann fluxes on only the fractures is not possible
264 // However, it is difficult to interpret this here as a fracture ending on the boundary
265 // could also be connected to a facet which is both boundary and fracture at the same time!
266 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
267 // we would have to find only those fractures that are at the boundary and aren't connected
268 // to a fracture which is a boundary.
269 LocalIndexType scvLocalIdx = element.subEntities(dim);
270 for (const auto& intersection : intersections(this->gridView(), element))
271 {
272 // first, obtain all vertex indices on this intersection
273 const auto& isGeometry = intersection.geometry();
274 const auto numCorners = isGeometry.corners();
275 const auto idxInInside = intersection.indexInInside();
276
277 std::vector<GridIndexType> isVertexIndices(numCorners);
278 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
279 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
280 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
281 dim);
282 // maybe add boundary scvf
283 if (intersection.boundary() && !intersection.neighbor())
284 {
285 numScvf_ += isGeometry.corners();
286 numBoundaryScvf_ += isGeometry.corners();
287
288 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
289 {
290 // find the scvs this scvf is belonging to
291 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
292 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
293 scvfs_[eIdx].emplace_back(geometryHelper,
294 intersection,
295 isGeometry,
296 isScvfLocalIdx,
297 scvfLocalIdx++,
298 std::move(localScvIndices));
299 }
300
301 // add all vertices on the intersection to the set of
302 // boundary vertices
303 const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
304 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
305 {
306 const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
307 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
308 boundaryDofIndices_[vIdxGlobal] = true;
309 }
310 }
311 // make sure we have no periodic boundaries
312 else if (intersection.boundary() && intersection.neighbor())
313 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
314
315 // maybe add fracture scvs & scvfs
316 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
317 {
318 for (auto vIdx : isVertexIndices)
319 fractureDofIndices_[vIdx] = true;
320
321 // add fracture scv for each vertex of intersection
322 numScv_ += numCorners;
323 const auto curNumScvs = scvs_[eIdx].size();
324 scvs_[eIdx].reserve(curNumScvs+numCorners);
325 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
326 scvs_[eIdx].emplace_back(geometryHelper,
327 intersection,
328 isGeometry,
329 vIdxLocal,
330 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
331 scvLocalIdx++,
332 idxInInside,
333 eIdx,
334 isVertexIndices[vIdxLocal]);
335
336 // add fracture scvf for each edge of the intersection in 3d
337 if (dim == 3)
338 {
339 const auto& faceRefElement = referenceElement(isGeometry);
340 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
341 {
342 // inside/outside scv indices in face local node numbering
343 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
344 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
345
346 // add offset to get the right scv indices
347 std::for_each( localScvIndices.begin(),
348 localScvIndices.end(),
349 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
350
351 // add scvf
352 numScvf_++;
353 scvfs_[eIdx].emplace_back(geometryHelper,
354 intersection,
355 isGeometry,
356 edgeIdx,
357 scvfLocalIdx++,
358 std::move(localScvIndices),
359 intersection.boundary());
360 }
361 }
362
363 // dim == 2, intersection is an edge, make 1 scvf
364 else
365 {
366 // inside/outside scv indices in face local node numbering
367 std::vector<LocalIndexType> localScvIndices({0, 1});
368
369 // add offset such that the fracture scvs above are addressed
370 std::for_each( localScvIndices.begin(),
371 localScvIndices.end(),
372 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
373
374 // add scvf
375 numScvf_++;
376 scvfs_[eIdx].emplace_back(geometryHelper,
377 intersection,
378 isGeometry,
379 /*idxOnIntersection*/0,
380 scvfLocalIdx++,
381 std::move(localScvIndices),
382 intersection.boundary());
383 }
384 }
385 }
386 }
387 }
388
389 const FeCache feCache_;
390
391 std::vector<std::vector<SubControlVolume>> scvs_;
392 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
393
394 // TODO do we need those?
395 std::size_t numScv_;
396 std::size_t numScvf_;
397 std::size_t numBoundaryScvf_;
398
399 // vertices on the boundary & fracture facets
400 std::vector<bool> boundaryDofIndices_;
401 std::vector<bool> fractureDofIndices_;
402};
403
411template<class Scalar, class GV, class Traits>
412class BoxDfmFVGridGeometry<Scalar, GV, false, Traits>
413: public BaseGridGeometry<GV, Traits>
414{
417 using GridIndexType = typename GV::IndexSet::IndexType;
418
419 static const int dim = GV::dimension;
420 static const int dimWorld = GV::dimensionworld;
421
422 using Element = typename GV::template Codim<0>::Entity;
423 using Intersection = typename GV::Intersection;
424 using CoordScalar = typename GV::ctype;
425
426public:
429 static constexpr DiscretizationMethod discMethod{};
430
432 using LocalView = typename Traits::template LocalView<ThisType, false>;
434 using SubControlVolume = typename Traits::SubControlVolume;
436 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
440 using DofMapper = typename Traits::VertexMapper;
442 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
444 using GridView = GV;
447
449 template< class FractureGridAdapter >
450 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
451 : ParentType(gridView)
452 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
453 {
454 update_(fractureGridAdapter);
455 }
456
459 const DofMapper& dofMapper() const
460 { return this->vertexMapper(); }
461
463 std::size_t numScv() const
464 { return numScv_; }
465
467 std::size_t numScvf() const
468 { return numScvf_; }
469
472 std::size_t numBoundaryScvf() const
473 { return numBoundaryScvf_; }
474
476 std::size_t numDofs() const
477 { return this->gridView().size(dim); }
478
480 template< class FractureGridAdapter >
481 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
482 {
483 ParentType::update(gridView);
484 updateFacetMapper_();
485 update_(fractureGridAdapter);
486 }
487
489 template< class FractureGridAdapter >
490 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
491 {
492 ParentType::update(std::move(gridView));
493 updateFacetMapper_();
494 update_(fractureGridAdapter);
495 }
496
498 const FeCache& feCache() const { return feCache_; }
500 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
502 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
504 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
505
507 bool isOnFracture(const Element& element, const Intersection& intersection) const
508 { return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
509
511 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
512 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
513
514private:
515
516 void updateFacetMapper_()
517 {
518 facetMapper_.update(this->gridView());
519 }
520
521 template< class FractureGridAdapter >
522 void update_(const FractureGridAdapter& fractureGridAdapter)
523 {
524 boundaryDofIndices_.assign(numDofs(), false);
525 fractureDofIndices_.assign(numDofs(), false);
526 facetOnFracture_.assign(this->gridView().size(1), false);
527
528 // save global data on the grid's scvs and scvfs
529 // TODO do we need those information?
530 numScv_ = 0;
531 numScvf_ = 0;
532 numBoundaryScvf_ = 0;
533 for (const auto& element : elements(this->gridView()))
534 {
535 numScv_ += element.subEntities(dim);
536 numScvf_ += element.subEntities(dim-1);
537
538 const auto elementGeometry = element.geometry();
539 const auto refElement = referenceElement(elementGeometry);
540
541 // store the sub control volume face indices on the domain boundary
542 for (const auto& intersection : intersections(this->gridView(), element))
543 {
544 // first, obtain all vertex indices on this intersection
545 const auto& isGeometry = intersection.geometry();
546 const auto numCorners = isGeometry.corners();
547 const auto idxInInside = intersection.indexInInside();
548
549 std::vector<GridIndexType> isVertexIndices(numCorners);
550 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
551 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
552 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
553 dim);
554
555 if (intersection.boundary() && !intersection.neighbor())
556 {
557 numScvf_ += numCorners;
558 numBoundaryScvf_ += numCorners;
559
560 // add all vertices on the intersection to the set of
561 // boundary vertices
562 const auto fIdx = intersection.indexInInside();
563 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
564 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
565 {
566 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
567 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
568 boundaryDofIndices_[vIdxGlobal] = true;
569 }
570 }
571
572 // make sure we have no periodic boundaries
573 else if (intersection.boundary() && intersection.neighbor())
574 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
575
576 // maybe add fracture scvs & scvfs
577 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
578 {
579 facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = true;
580 for (auto vIdx : isVertexIndices)
581 fractureDofIndices_[vIdx] = true;
582
583 const auto isGeometry = intersection.geometry();
584 numScv_ += isGeometry.corners();
585 numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
586 }
587 }
588 }
589 }
590
591 const FeCache feCache_;
592
593 // Information on the global number of geometries
594 // TODO do we need those information?
595 std::size_t numScv_;
596 std::size_t numScvf_;
597 std::size_t numBoundaryScvf_;
598
599 // vertices on the boundary & fracture facets
600 std::vector<bool> boundaryDofIndices_;
601 std::vector<bool> fractureDofIndices_;
602
603 // facet mapper and markers which facets lie on interior boundaries
604 typename Traits::FacetMapper facetMapper_;
605 std::vector<bool> facetOnFracture_;
606};
607
608} // end namespace Dumux
609
610#endif
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:41
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:414
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:438
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:500
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:504
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:450
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:463
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:507
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:476
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:472
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:481
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:440
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:432
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:467
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:442
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:498
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:459
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:436
Detail::BoxDfmGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:446
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:434
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:502
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:490
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:511
Base class for the finite volume geometry vector for box schemes that consider extra connectivity bet...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:98
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:117
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:125
Detail::BoxDfmGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:129
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:182
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:158
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:123
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:133
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:171
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:180
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:121
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:115
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:141
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:149
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:145
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:191
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:154
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:184
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:178
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:163
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:119
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:188
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:186
Base class for the finite volume geometry vector for box schemes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:83
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:47
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:65
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
Defines the default element and vertex mapper types.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Dune::Std::detected_or_t< Dumux::BoxDfmGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxDfmGeometryHelper_t
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:46
typename T::GeometryHelper SpecifiesGeometryHelper
Definition: basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: common/pdesolver.hh:24
Base class for the local finite volume geometry for the box discrete fracture model.
Helper class constructing the dual grid finite volume geometries 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.
The default traits for the box finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:60
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:68
Definition: method.hh:46