version 3.11-dev
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
33
34#include "fvelementgeometry.hh"
35#include "geometryhelper.hh"
36#include "subcontrolvolume.hh"
38
39namespace Dumux {
40
41namespace Detail {
42template<class GV, class T>
43using BoxDfmGeometryHelper_t = Dune::Std::detected_or_t<
46 T
47>;
48} // end namespace Detail
49
58template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
60: public MapperTraits
61{
64
65 template<class GridGeometry, bool enableCache>
67
68 // Mapper type for mapping edges
69 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
70};
71
80template<class Scalar,
81 class GridView,
82 bool enableGridGeometryCache = false,
85
96template<class Scalar, class GV, class Traits>
97class BoxDfmFVGridGeometry<Scalar, GV, true, Traits>
98: public BaseGridGeometry<GV, Traits>
99{
102 using GridIndexType = typename IndexTraits<GV>::GridIndex;
103 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
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
111public:
114 static constexpr DiscretizationMethod discMethod{};
115
117 using LocalView = typename Traits::template LocalView<ThisType, true>;
119 using SubControlVolume = typename Traits::SubControlVolume;
121 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
125 using DofMapper = typename Traits::VertexMapper;
127 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
129 using GridView = GV;
130
132 template< class FractureGridAdapter >
133 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
134 : ParentType(gridView)
135 , cache_(*this)
136 {
137 update_(fractureGridAdapter);
138 }
139
142 const DofMapper& dofMapper() const
143 { return this->vertexMapper(); }
144
146 std::size_t numScv() const
147 { return numScv_; }
148
150 std::size_t numScvf() const
151 { return numScvf_; }
152
155 std::size_t numBoundaryScvf() const
156 { return numBoundaryScvf_; }
157
159 std::size_t numDofs() const
160 { return this->gridView().size(dim); }
161
163 template< class FractureGridAdapter >
164 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
165 {
166 ParentType::update(gridView);
167 update_(fractureGridAdapter);
168 }
169
171 template< class FractureGridAdapter >
172 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
173 {
174 ParentType::update(std::move(gridView));
175 update_(fractureGridAdapter);
176 }
177
179 const FeCache& feCache() const { return feCache_; }
181 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
183 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
185 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
186
188 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
189 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
190
192 friend inline LocalView localView(const BoxDfmFVGridGeometry& gg)
193 { return { gg.cache_ }; }
194
195private:
196
197 class BoxDfmGridGeometryCache
198 {
199 friend class BoxDfmFVGridGeometry;
200 public:
202 using GeometryHelper = Detail::BoxDfmGeometryHelper_t<GV, Traits>;
203
204 explicit BoxDfmGridGeometryCache(const BoxDfmFVGridGeometry& gg)
205 : gridGeometry_(&gg)
206 {}
207
208 const BoxDfmFVGridGeometry& gridGeometry() const
209 { return *gridGeometry_; }
210
212 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
213 { return scvs_[eIdx]; }
214
216 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
217 { return scvfs_[eIdx]; }
218
219
220 private:
221 void clear_()
222 {
223 scvs_.clear();
224 scvfs_.clear();
225 }
226
227 std::vector<std::vector<SubControlVolume>> scvs_;
228 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
229
230 const BoxDfmFVGridGeometry* gridGeometry_;
231 };
232
233public:
236 using Cache = BoxDfmGridGeometryCache;
237
238private:
239 using GeometryHelper = typename Cache::GeometryHelper;
240
241 template< class FractureGridAdapter >
242 void update_(const FractureGridAdapter& fractureGridAdapter)
243 {
244 cache_.clear_();
245
246 const auto numElements = this->gridView().size(0);
247 cache_.scvs_.resize(numElements);
248 cache_.scvfs_.resize(numElements);
249
250 boundaryDofIndices_.assign(numDofs(), false);
251 fractureDofIndices_.assign(this->gridView.size(dim), false);
252
253 numScv_ = 0;
254 numScvf_ = 0;
255 numBoundaryScvf_ = 0;
256 // Build the SCV and SCV faces
257 for (const auto& element : elements(this->gridView()))
258 {
259 // fill the element map with seeds
260 const auto eIdx = this->elementMapper().index(element);
261
262 // count
263 numScv_ += element.subEntities(dim);
264 numScvf_ += element.subEntities(dim-1);
265
266 // get the element geometry
267 auto elementGeometry = element.geometry();
268 const auto refElement = referenceElement(elementGeometry);
269
270 // instantiate the geometry helper
271 GeometryHelper geometryHelper(elementGeometry);
272
273 // construct the sub control volumes
274 cache_.scvs_[eIdx].resize(elementGeometry.corners());
275 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
276 {
277 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
278
279 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
280 scvLocalIdx,
281 eIdx,
282 dofIdxGlobal);
283 }
284
285 // construct the sub control volume faces
286 LocalIndexType scvfLocalIdx = 0;
287 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
288 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
289 {
290 // find the global and local scv indices this scvf is belonging to
291 std::array<LocalIndexType, 2> localScvIndices{{
292 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
293 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
294 }};
295
296 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
297 element,
298 scvfLocalIdx,
299 std::move(localScvIndices));
300 }
301
302 // construct the ...
303 // ... sub-control volume faces on the domain boundary
304 // ... sub-control volumes on fracture facets
305 // ... sub-control volume faces on fracture facets
306 // NOTE We do not construct fracture scvfs on boundaries here!
307 // That means specifying Neumann fluxes on only the fractures is not possible
308 // However, it is difficult to interpret this here as a fracture ending on the boundary
309 // could also be connected to a facet which is both boundary and fracture at the same time!
310 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
311 // we would have to find only those fractures that are at the boundary and aren't connected
312 // to a fracture which is a boundary.
313 LocalIndexType scvLocalIdx = element.subEntities(dim);
314 for (const auto& intersection : intersections(this->gridView(), element))
315 {
316 // first, obtain all vertex indices on this intersection
317 const auto& isGeometry = intersection.geometry();
318 const auto numCorners = isGeometry.corners();
319 const auto idxInInside = intersection.indexInInside();
320
321 std::vector<GridIndexType> isVertexIndices(numCorners);
322 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
323 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
324 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
325 dim);
326 // maybe add boundary scvf
327 if (intersection.boundary() && !intersection.neighbor())
328 {
329 numScvf_ += isGeometry.corners();
330 numBoundaryScvf_ += isGeometry.corners();
331
332 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
333 {
334 // find the scvs this scvf is belonging to
335 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
336 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
337 cache_.scvfs_[eIdx].emplace_back(geometryHelper,
338 intersection,
339 isScvfLocalIdx,
340 scvfLocalIdx++,
341 std::move(localScvIndices));
342 }
343
344 // add all vertices on the intersection to the set of
345 // boundary vertices
346 const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
347 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
348 {
349 const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
350 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
351 boundaryDofIndices_[vIdxGlobal] = true;
352 }
353 }
354 // make sure we have no periodic boundaries
355 else if (intersection.boundary() && intersection.neighbor())
356 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
357
358 // maybe add fracture scvs & scvfs
359 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
360 {
361 for (auto vIdx : isVertexIndices)
362 fractureDofIndices_[vIdx] = true;
363
364 // add fracture scv for each vertex of intersection
365 numScv_ += numCorners;
366 const auto curNumScvs = cache_.scvs_[eIdx].size();
367 cache_.scvs_[eIdx].reserve(curNumScvs+numCorners);
368 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
369 cache_.scvs_[eIdx].emplace_back(geometryHelper,
370 intersection,
371 isGeometry,
372 vIdxLocal,
373 static_cast<LocalIndexType>(refElement.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& faceRefElement = referenceElement(isGeometry);
383 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
384 {
385 // inside/outside scv indices in face local node numbering
386 std::array<LocalIndexType, 2> localScvIndices{{
387 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
388 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))
389 }};
390
391 // add offset to get the right scv indices
392 std::for_each( localScvIndices.begin(),
393 localScvIndices.end(),
394 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
395
396 // add scvf
397 numScvf_++;
398 cache_.scvfs_[eIdx].emplace_back(geometryHelper,
399 intersection,
400 edgeIdx,
401 scvfLocalIdx++,
402 std::move(localScvIndices),
403 intersection.boundary());
404 }
405 }
406
407 // dim == 2, intersection is an edge, make 1 scvf
408 else
409 {
410 // inside/outside scv indices in face local node numbering
411 std::array<LocalIndexType, 2> localScvIndices{{0, 1}};
412
413 // add offset such that the fracture scvs above are addressed
414 std::for_each( localScvIndices.begin(),
415 localScvIndices.end(),
416 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
417
418 // add scvf
419 numScvf_++;
420 cache_.scvfs_[eIdx].emplace_back(geometryHelper,
421 intersection,
422 /*idxOnIntersection*/0,
423 scvfLocalIdx++,
424 std::move(localScvIndices),
425 intersection.boundary());
426 }
427 }
428 }
429 }
430 }
431
432 const FeCache feCache_;
433
434 // TODO do we need those?
435 std::size_t numScv_;
436 std::size_t numScvf_;
437 std::size_t numBoundaryScvf_;
438
439 // vertices on the boundary & fracture facets
440 std::vector<bool> boundaryDofIndices_;
441 std::vector<bool> fractureDofIndices_;
442
443 Cache cache_;
444};
445
453template<class Scalar, class GV, class Traits>
454class BoxDfmFVGridGeometry<Scalar, GV, false, Traits>
455: public BaseGridGeometry<GV, Traits>
456{
459 using GridIndexType = typename GV::IndexSet::IndexType;
460
461 static const int dim = GV::dimension;
462 static const int dimWorld = GV::dimensionworld;
463
464 using Element = typename GV::template Codim<0>::Entity;
465 using Intersection = typename GV::Intersection;
466 using CoordScalar = typename GV::ctype;
467
468public:
471 static constexpr DiscretizationMethod discMethod{};
472
474 using LocalView = typename Traits::template LocalView<ThisType, false>;
476 using SubControlVolume = typename Traits::SubControlVolume;
478 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
482 using DofMapper = typename Traits::VertexMapper;
484 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
486 using GridView = GV;
487
489 template< class FractureGridAdapter >
490 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
491 : ParentType(gridView)
492 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
493 , cache_(*this)
494 {
495 update_(fractureGridAdapter);
496 }
497
500 const DofMapper& dofMapper() const
501 { return this->vertexMapper(); }
502
504 std::size_t numScv() const
505 { return numScv_; }
506
508 std::size_t numScvf() const
509 { return numScvf_; }
510
513 std::size_t numBoundaryScvf() const
514 { return numBoundaryScvf_; }
515
517 std::size_t numDofs() const
518 { return this->gridView().size(dim); }
519
521 template< class FractureGridAdapter >
522 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
523 {
524 ParentType::update(gridView);
525 updateFacetMapper_();
526 update_(fractureGridAdapter);
527 }
528
530 template< class FractureGridAdapter >
531 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
532 {
533 ParentType::update(std::move(gridView));
534 updateFacetMapper_();
535 update_(fractureGridAdapter);
536 }
537
539 const FeCache& feCache() const { return feCache_; }
541 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
543 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
545 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
546
548 bool isOnFracture(const Element& element, const Intersection& intersection) const
549 { return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
550
552 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
553 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
554
556 friend inline LocalView localView(const BoxDfmFVGridGeometry& gg)
557 { return { gg.cache_ }; }
558
559private:
560 class BoxDfmGridGeometryCache
561 {
562 friend class BoxDfmFVGridGeometry;
563 public:
565 using GeometryHelper = Detail::BoxDfmGeometryHelper_t<GV, Traits>;
566
567 explicit BoxDfmGridGeometryCache(const BoxDfmFVGridGeometry& gg)
568 : gridGeometry_(&gg)
569 {}
570
571 const BoxDfmFVGridGeometry& gridGeometry() const
572 { return *gridGeometry_; }
573
574 private:
575 const BoxDfmFVGridGeometry* gridGeometry_;
576 };
577
578public:
581 using Cache = BoxDfmGridGeometryCache;
582
583private:
584
585 void updateFacetMapper_()
586 {
587 facetMapper_.update(this->gridView());
588 }
589
590 template< class FractureGridAdapter >
591 void update_(const FractureGridAdapter& fractureGridAdapter)
592 {
593 boundaryDofIndices_.assign(numDofs(), false);
594 fractureDofIndices_.assign(numDofs(), false);
595 facetOnFracture_.assign(this->gridView().size(1), false);
596
597 // save global data on the grid's scvs and scvfs
598 // TODO do we need those information?
599 numScv_ = 0;
600 numScvf_ = 0;
601 numBoundaryScvf_ = 0;
602 for (const auto& element : elements(this->gridView()))
603 {
604 numScv_ += element.subEntities(dim);
605 numScvf_ += element.subEntities(dim-1);
606
607 const auto elementGeometry = element.geometry();
608 const auto refElement = referenceElement(elementGeometry);
609
610 // store the sub control volume face indices on the domain boundary
611 for (const auto& intersection : intersections(this->gridView(), element))
612 {
613 // first, obtain all vertex indices on this intersection
614 const auto& isGeometry = intersection.geometry();
615 const auto numCorners = isGeometry.corners();
616 const auto idxInInside = intersection.indexInInside();
617
618 std::vector<GridIndexType> isVertexIndices(numCorners);
619 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
620 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
621 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
622 dim);
623
624 if (intersection.boundary() && !intersection.neighbor())
625 {
626 numScvf_ += numCorners;
627 numBoundaryScvf_ += numCorners;
628
629 // add all vertices on the intersection to the set of
630 // boundary vertices
631 const auto fIdx = intersection.indexInInside();
632 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
633 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
634 {
635 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
636 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
637 boundaryDofIndices_[vIdxGlobal] = true;
638 }
639 }
640
641 // make sure we have no periodic boundaries
642 else if (intersection.boundary() && intersection.neighbor())
643 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
644
645 // maybe add fracture scvs & scvfs
646 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
647 {
648 facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = true;
649 for (auto vIdx : isVertexIndices)
650 fractureDofIndices_[vIdx] = true;
651
652 const auto isGeometry = intersection.geometry();
653 numScv_ += isGeometry.corners();
654 numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
655 }
656 }
657 }
658 }
659
660 const FeCache feCache_;
661
662 // Information on the global number of geometries
663 // TODO do we need those information?
664 std::size_t numScv_;
665 std::size_t numScvf_;
666 std::size_t numBoundaryScvf_;
667
668 // vertices on the boundary & fracture facets
669 std::vector<bool> boundaryDofIndices_;
670 std::vector<bool> fractureDofIndices_;
671
672 // facet mapper and markers which facets lie on interior boundaries
673 typename Traits::FacetMapper facetMapper_;
674 std::vector<bool> facetOnFracture_;
675
676 Cache cache_;
677};
678
679} // end namespace Dumux
680
681#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:44
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:456
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:480
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:541
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:545
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:490
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:504
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:548
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:517
BoxDfmGridGeometryCache Cache
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:581
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:513
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:522
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:482
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:474
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:508
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:484
friend LocalView localView(const BoxDfmFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:556
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:539
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:500
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:478
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:476
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:543
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:531
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:552
Base class for the finite volume geometry vector for box schemes that consider extra connectivity bet...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:99
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:119
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:127
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:159
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:125
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:133
BoxDfmGridGeometryCache Cache
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:236
friend LocalView localView(const BoxDfmFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:192
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:172
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:123
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:117
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:142
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:150
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:146
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:188
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:155
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:181
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:179
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:164
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:121
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:185
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:183
Base class for the finite volume geometry vector for box schemes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:84
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:63
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:62
Defines the default element and vertex mapper types.
Helper classes to compute the integration elements.
Defines the index types used for grid and local indices.
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:47
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
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
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:61
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:69
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26