version 3.11-dev
Loading...
Searching...
No Matches
multidomain/facet/box/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//
15#ifndef DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
16#define DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
17
18#include <algorithm>
19#include <utility>
20
21#include <dune/grid/common/mcmgmapper.hh>
22#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
23
30
34
35namespace Dumux {
36
37namespace Detail {
38template<class GV, class T>
39using BoxFacetCouplingGeometryHelper_t = Dune::Std::detected_or_t<
42 T
43>;
44} // end namespace Detail
45
53template<class GridView>
55{
56 // use a specialized version of the box scvf
59
60 template<class GridGeometry, bool enableCache>
62
63 // per default we use an mcmg mapper for the elements
64 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
65 // the default vertex mapper is the enriched vertex dof mapper
67 // Mapper type for mapping edges
68 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
69};
70
78template<class Scalar,
79 class GridView,
80 bool enableGridGeometryCache = false,
83
91template<class Scalar, class GV, class Traits>
92class BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits>
93: public BaseGridGeometry<GV, Traits>
94{
96 using ParentType = BaseGridGeometry<GV, Traits>;
97 using GridIndexType = typename IndexTraits<GV>::GridIndex;
98 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
99
100 using Element = typename GV::template Codim<0>::Entity;
101 using CoordScalar = typename GV::ctype;
102 static const int dim = GV::dimension;
103 static const int dimWorld = GV::dimensionworld;
104
105public:
109
111 using LocalView = typename Traits::template LocalView<ThisType, true>;
113 using SubControlVolume = typename Traits::SubControlVolume;
115 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
119 using DofMapper = typename Traits::VertexMapper;
121 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
123 using GridView = GV;
124
126 template<class FacetGridView, class CodimOneGridAdapter>
128 const FacetGridView& facetGridView,
129 const CodimOneGridAdapter& codimOneGridAdapter,
130 bool verbose = false)
131 : ParentType(gridView)
132 , cache_(*this)
133 {
134 update_(facetGridView, codimOneGridAdapter, verbose);
135 }
136
138 const DofMapper& dofMapper() const
139 { return this->vertexMapper(); }
140
142 std::size_t numScv() const
143 { return numScv_; }
144
146 std::size_t numScvf() const
147 { return numScvf_; }
148
151 std::size_t numBoundaryScvf() const
152 { return numBoundaryScvf_; }
153
155 std::size_t numDofs() const
156 { return this->vertexMapper().size(); }
157
171 template<class FacetGridView, class CodimOneGridAdapter>
173 const FacetGridView& facetGridView,
174 const CodimOneGridAdapter& codimOneGridAdapter,
175 bool verbose = false)
176 {
178 update_(facetGridView, codimOneGridAdapter, verbose);
179 }
180
182 template<class FacetGridView, class CodimOneGridAdapter>
184 const FacetGridView& facetGridView,
185 const CodimOneGridAdapter& codimOneGridAdapter,
186 bool verbose = false)
187 {
188 ParentType::update(std::move(gridView));
189 update_(facetGridView, codimOneGridAdapter, verbose);
190 }
191
193 const FeCache& feCache() const
194 { return feCache_; }
195
197 bool dofOnBoundary(GridIndexType dofIdx) const
198 { return boundaryDofIndices_[dofIdx]; }
199
201 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
202 { return interiorBoundaryDofIndices_[dofIdx]; }
203
205 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
206 { return false; }
207
209 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
210 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
211
214 { return { gg.cache_ }; }
215
216private:
217 class BoxFacetCouplingGridGeometryCache
218 {
220 public:
223
224 explicit BoxFacetCouplingGridGeometryCache(const BoxFacetCouplingFVGridGeometry& gg)
225 : gridGeometry_(&gg)
226 {}
227
228 const BoxFacetCouplingFVGridGeometry& gridGeometry() const
229 { return *gridGeometry_; }
230
232 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
233 { return scvs_[eIdx]; }
234
236 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
237 { return scvfs_[eIdx]; }
238
239
240 private:
241 void clear_()
242 {
243 scvs_.clear();
244 scvfs_.clear();
245 }
246
247 std::vector<std::vector<SubControlVolume>> scvs_;
248 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
249
250 const BoxFacetCouplingFVGridGeometry* gridGeometry_;
251 };
252
253public:
256 using Cache = BoxFacetCouplingGridGeometryCache;
257
258private:
259 using GeometryHelper = typename Cache::GeometryHelper;
260
261 template<class FacetGridView, class CodimOneGridAdapter>
262 void update_(const FacetGridView& facetGridView,
263 const CodimOneGridAdapter& codimOneGridAdapter,
264 bool verbose = false)
265 {
266 // enrich the vertex mapper subject to the provided facet grid
267 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
268
269 cache_.clear_();
270
271 // resize containers
272 const auto numDof = numDofs();
273 const auto numElements = this->gridView().size(0);
274 cache_.scvs_.resize(numElements);
275 cache_.scvfs_.resize(numElements);
276 boundaryDofIndices_.assign(numDof, false);
277 interiorBoundaryDofIndices_.assign(numDof, false);
278
279 // Build the SCV and SCV faces
280 numScv_ = 0;
281 numScvf_ = 0;
282 numBoundaryScvf_ = 0;
283 for (const auto& element : elements(this->gridView()))
284 {
285 auto eIdx = this->elementMapper().index(element);
286
287 // keep track of number of scvs and scvfs
288 numScv_ += element.subEntities(dim);
289 numScvf_ += element.subEntities(dim-1);
290
291 // get the element geometry
292 auto elementGeometry = element.geometry();
293 const auto refElement = referenceElement(elementGeometry);
294
295 // instantiate the geometry helper
296 GeometryHelper geometryHelper(elementGeometry);
297
298 // construct the sub control volumes
299 cache_.scvs_[eIdx].resize(elementGeometry.corners());
300 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
301 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper.getScvCorners(scvLocalIdx),
302 scvLocalIdx,
303 eIdx,
304 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
305
306 // construct the sub control volume faces
307 LocalIndexType scvfLocalIdx = 0;
308 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
309 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
310 {
311 // find the global and local scv indices this scvf is belonging to
312 std::array<LocalIndexType, 2> localScvIndices{{
313 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
314 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
315 }};
316
317 // create the sub-control volume face
318 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
319 element,
320 scvfLocalIdx,
321 std::move(localScvIndices));
322 }
323
324 // construct the sub control volume faces on the domain/interior boundaries
325 // skip handled facets (necessary for e.g. Dune::FoamGrid)
326 std::vector<unsigned int> handledFacets;
327 for (const auto& intersection : intersections(this->gridView(), element))
328 {
329 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
330 continue;
331
332 handledFacets.push_back(intersection.indexInInside());
333
334 // determine if all corners live on the facet grid
335 const auto isGeometry = intersection.geometry();
336 const auto numFaceCorners = isGeometry.corners();
337 const auto idxInInside = intersection.indexInInside();
338 const auto boundary = intersection.boundary();
339
340 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
341 for (int i = 0; i < numFaceCorners; ++i)
342 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
343
344 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
345 for (int i = 0; i < numFaceCorners; ++i)
346 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
347
348 // if the vertices compose a facet element, the intersection is on facet grid
349 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
350
351 // make sure there are no periodic boundaries
352 if (boundary && intersection.neighbor())
353 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
354
355 if (isOnFacet || boundary)
356 {
357 // keep track of number of faces
358 numScvf_ += numFaceCorners;
359 numBoundaryScvf_ += int(boundary)*numFaceCorners;
360
361 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
362 {
363 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
364 std::array<LocalIndexType, 2> localScvIndices{{
365 vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]
366 }};
367
368 // create the sub-control volume face
369 cache_.scvfs_[eIdx].emplace_back(geometryHelper,
370 intersection,
371 isScvfLocalIdx,
372 scvfLocalIdx,
373 std::move(localScvIndices),
374 boundary,
375 isOnFacet);
376
377 // Mark vertices to be on domain and/or interior boundary
378 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
379 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
380 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
381
382 // increment local counter
383 scvfLocalIdx++;
384 }
385 }
386 }
387 }
388 }
389
390 const FeCache feCache_;
391
392 // TODO do we need those?
393 std::size_t numScv_;
394 std::size_t numScvf_;
395 std::size_t numBoundaryScvf_;
396
397 // vertices on domain/interior boundaries
398 std::vector<bool> boundaryDofIndices_;
399 std::vector<bool> interiorBoundaryDofIndices_;
400
401 Cache cache_;
402};
403
411template<class Scalar, class GV, class Traits>
412class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
413: public BaseGridGeometry<GV, Traits>
414{
416 using ParentType = BaseGridGeometry<GV, Traits>;
417 using GridIndexType = typename IndexTraits<GV>::GridIndex;
418 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
419
420 static const int dim = GV::dimension;
421 static const int dimWorld = GV::dimensionworld;
422
423 using Element = typename GV::template Codim<0>::Entity;
424 using Intersection = typename GV::Intersection;
425 using CoordScalar = typename GV::ctype;
426
427public:
431
433 using LocalView = typename Traits::template LocalView<ThisType, false>;
435 using SubControlVolume = typename Traits::SubControlVolume;
437 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
441 using DofMapper = typename Traits::VertexMapper;
443 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
445 using GridView = GV;
446
448 template<class FacetGridView, class CodimOneGridAdapter>
450 const FacetGridView& facetGridView,
451 const CodimOneGridAdapter& codimOneGridAdapter,
452 bool verbose = false)
453 : ParentType(gridView)
454 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
455 , cache_(*this)
456 {
457 update_(facetGridView, codimOneGridAdapter, verbose);
458 }
459
462 const DofMapper& dofMapper() const
463 { return this->vertexMapper(); }
464
466 std::size_t numScv() const
467 { return numScv_; }
468
470 std::size_t numScvf() const
471 { return numScvf_; }
472
475 std::size_t numBoundaryScvf() const
476 { return numBoundaryScvf_; }
477
479 std::size_t numDofs() const
480 { return this->vertexMapper().size(); }
481
495 template<class FacetGridView, class CodimOneGridAdapter>
497 const FacetGridView& facetGridView,
498 const CodimOneGridAdapter& codimOneGridAdapter,
499 bool verbose = false)
500 {
502 updateFacetMapper_();
503 update_(facetGridView, codimOneGridAdapter, verbose);
504 }
505
507 template<class FacetGridView, class CodimOneGridAdapter>
509 const FacetGridView& facetGridView,
510 const CodimOneGridAdapter& codimOneGridAdapter,
511 bool verbose = false)
512 {
513 ParentType::update(std::move(gridView));
514 updateFacetMapper_();
515 update_(facetGridView, codimOneGridAdapter, verbose);
516 }
517
519 const FeCache& feCache() const
520 { return feCache_; }
521
523 bool dofOnBoundary(unsigned int dofIdx) const
524 { return boundaryDofIndices_[dofIdx]; }
525
527 bool dofOnInteriorBoundary(unsigned int dofIdx) const
528 { return interiorBoundaryDofIndices_[dofIdx]; }
529
531 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
532 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
533
535 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
536 { return false; }
537
539 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
540 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
541
544 { return { gg.cache_ }; }
545
546private:
547
548 class BoxFacetCouplingGridGeometryCache
549 {
551 public:
554
555 explicit BoxFacetCouplingGridGeometryCache(const BoxFacetCouplingFVGridGeometry& gg)
556 : gridGeometry_(&gg)
557 {}
558
559 const BoxFacetCouplingFVGridGeometry& gridGeometry() const
560 { return *gridGeometry_; }
561
562 private:
563 const BoxFacetCouplingFVGridGeometry* gridGeometry_;
564 };
565
566public:
569 using Cache = BoxFacetCouplingGridGeometryCache;
570
571private:
572
573 void updateFacetMapper_()
574 {
575 facetMapper_.update(this->gridView());
576 }
577
578 template<class FacetGridView, class CodimOneGridAdapter>
579 void update_(const FacetGridView& facetGridView,
580 const CodimOneGridAdapter& codimOneGridAdapter,
581 bool verbose)
582 {
583 // enrich the vertex mapper subject to the provided facet grid
584 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
585
586 // save global data on the grid's scvs and scvfs
587 // TODO do we need those information?
588 numScv_ = 0;
589 numScvf_ = 0;
590 numBoundaryScvf_ = 0;
591
592 const auto numDof = numDofs();
593 boundaryDofIndices_.assign(numDof, false);
594 interiorBoundaryDofIndices_.assign(numDof, false);
595 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
596 for (const auto& element : elements(this->gridView()))
597 {
598 numScv_ += element.subEntities(dim);
599 numScvf_ += element.subEntities(dim-1);
600
601 const auto elementGeometry = element.geometry();
602 const auto refElement = referenceElement(elementGeometry);
603
604 // store the sub control volume face indices on the domain/interior boundary
605 // skip handled facets (necessary for e.g. Dune::FoamGrid)
606 std::vector<unsigned int> handledFacets;
607 for (const auto& intersection : intersections(this->gridView(), element))
608 {
609 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
610 continue;
611
612 handledFacets.push_back(intersection.indexInInside());
613
614 // determine if all corners live on the facet grid
615 const auto isGeometry = intersection.geometry();
616 const auto numFaceCorners = isGeometry.corners();
617 const auto idxInInside = intersection.indexInInside();
618 const auto boundary = intersection.boundary();
619
620 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
621 for (int i = 0; i < numFaceCorners; ++i)
622 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
623
624 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
625 for (int i = 0; i < numFaceCorners; ++i)
626 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
627
628 // if all vertices are living on the facet grid, this is an interiour boundary
629 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
630
631 // make sure there are no periodic boundaries
632 if (boundary && intersection.neighbor())
633 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
634
635 if (isOnFacet || boundary)
636 {
637 numScvf_ += numFaceCorners;
638 numBoundaryScvf_ += int(boundary)*numFaceCorners;
639
640 // Mark vertices to be on domain and/or interior boundary
641 for (int i = 0; i < numFaceCorners; ++i)
642 {
643 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
644 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
645 if (isOnFacet)
646 {
647 interiorBoundaryDofIndices_[ dofIndex ] = true;
648 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
649 }
650 }
651 }
652 }
653 }
654 }
655
656 const FeCache feCache_;
657
658 // Information on the global number of geometries
659 // TODO do we need those information?
660 std::size_t numScv_;
661 std::size_t numScvf_;
662 std::size_t numBoundaryScvf_;
663
664 // vertices on domain/interior boundaries
665 std::vector<bool> boundaryDofIndices_;
666 std::vector<bool> interiorBoundaryDofIndices_;
667
668 // facet mapper and markers which facets lie on interior boundaries
669 typename Traits::FacetMapper facetMapper_;
670 std::vector<bool> facetIsOnInteriorBoundary_;
671
672 Cache cache_;
673};
674
675} // end namespace Dumux
676
677#endif
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices for constant grids.
Definition basegridgeometry.hh:112
Element element(GridIndexType eIdx) const
Get an element from a global element index.
Definition basegridgeometry.hh:142
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices for constant grids.
Definition basegridgeometry.hh:106
const GridView & gridView() const
Return the gridView this grid geometry object lives on.
Definition basegridgeometry.hh:100
void update(const GridView &gridView)
Update all fvElementGeometries (call this after grid adaption).
Definition basegridgeometry.hh:88
Base class for the element-local finite volume geometry for box models in the context of models consi...
Definition multidomain/facet/box/fvelementgeometry.hh:36
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition multidomain/facet/box/fvgridgeometry.hh:435
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition multidomain/facet/box/fvgridgeometry.hh:449
std::size_t numScv() const
The total number of sub control volumes.
Definition multidomain/facet/box/fvgridgeometry.hh:466
std::size_t numScvf() const
The total number of sun control volume faces.
Definition multidomain/facet/box/fvgridgeometry.hh:470
GV GridView
export the grid view type
Definition multidomain/facet/box/fvgridgeometry.hh:445
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition multidomain/facet/box/fvgridgeometry.hh:429
std::size_t numDofs() const
The total number of degrees of freedom.
Definition multidomain/facet/box/fvgridgeometry.hh:479
void update(GridView &&gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition multidomain/facet/box/fvgridgeometry.hh:508
std::size_t numBoundaryScvf() const
Definition multidomain/facet/box/fvgridgeometry.hh:475
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition multidomain/facet/box/fvgridgeometry.hh:443
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition multidomain/facet/box/fvgridgeometry.hh:519
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition multidomain/facet/box/fvgridgeometry.hh:441
const DofMapper & dofMapper() const
Definition multidomain/facet/box/fvgridgeometry.hh:462
friend LocalView localView(const BoxFacetCouplingFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition multidomain/facet/box/fvgridgeometry.hh:543
static constexpr DiscretizationMethod discMethod
Definition multidomain/facet/box/fvgridgeometry.hh:430
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition multidomain/facet/box/fvgridgeometry.hh:433
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:539
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition multidomain/facet/box/fvgridgeometry.hh:535
BoxFacetCouplingGridGeometryCache Cache
Definition multidomain/facet/box/fvgridgeometry.hh:569
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition multidomain/facet/box/fvgridgeometry.hh:437
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition multidomain/facet/box/fvgridgeometry.hh:439
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:523
void update(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition multidomain/facet/box/fvgridgeometry.hh:496
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:527
bool isOnInteriorBoundary(const Element &element, const Intersection &intersection) const
returns true if an intersection is on an interior boundary
Definition multidomain/facet/box/fvgridgeometry.hh:531
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:201
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:209
void update(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition multidomain/facet/box/fvgridgeometry.hh:172
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition multidomain/facet/box/fvgridgeometry.hh:121
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition multidomain/facet/box/fvgridgeometry.hh:117
static constexpr DiscretizationMethod discMethod
Definition multidomain/facet/box/fvgridgeometry.hh:108
BoxFacetCouplingGridGeometryCache Cache
Definition multidomain/facet/box/fvgridgeometry.hh:256
GV GridView
export the grid view type
Definition multidomain/facet/box/fvgridgeometry.hh:123
std::size_t numScv() const
The total number of sub control volumes.
Definition multidomain/facet/box/fvgridgeometry.hh:142
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition multidomain/facet/box/fvgridgeometry.hh:205
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition multidomain/facet/box/fvgridgeometry.hh:197
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition multidomain/facet/box/fvgridgeometry.hh:119
std::size_t numBoundaryScvf() const
Definition multidomain/facet/box/fvgridgeometry.hh:151
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition multidomain/facet/box/fvgridgeometry.hh:138
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition multidomain/facet/box/fvgridgeometry.hh:113
friend LocalView localView(const BoxFacetCouplingFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition multidomain/facet/box/fvgridgeometry.hh:213
std::size_t numScvf() const
The total number of sun control volume faces.
Definition multidomain/facet/box/fvgridgeometry.hh:146
void update(GridView &&gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (call this after grid adaption)
Definition multidomain/facet/box/fvgridgeometry.hh:183
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition multidomain/facet/box/fvgridgeometry.hh:193
std::size_t numDofs() const
The total number of degrees of freedom.
Definition multidomain/facet/box/fvgridgeometry.hh:155
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition multidomain/facet/box/fvgridgeometry.hh:127
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition multidomain/facet/box/fvgridgeometry.hh:111
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition multidomain/facet/box/fvgridgeometry.hh:107
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition multidomain/facet/box/fvgridgeometry.hh:115
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition multidomain/facet/box/fvgridgeometry.hh:82
Class for a sub control volume face in the box method, i.e a part of the boundary of a sub control vo...
Definition multidomain/facet/box/subcontrolvolumeface.hh:39
Create sub control volumes and sub control volume face geometries.
Definition boxgeometryhelper.hh:257
the sub control volume for the box scheme
Definition discretization/box/subcontrolvolume.hh:58
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition codimonegridadapter.hh:40
bool composeFacetElement(const IndexStorage &bulkVertexIndices) const
Returns true if a given set of bulk vertex indices make up a facet grid element.
Definition codimonegridadapter.hh:180
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...
Definition vertexmapper.hh:126
the sub control volume for the box scheme
Helper classes to compute the integration elements.
BaseGridGeometry(std::shared_ptr< BaseImplementation > impl)
Constructor from a BaseImplementation.
Definition basegridgeometry.hh:72
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Base class for the element-local finite volume geometry for box models in the context of models consi...
Base class for a sub control volume face of the box method in the context of of models considering co...
Definition cvfelocalresidual.hh:25
Dune::Std::detected_or_t< Dumux::BoxGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxFacetCouplingGeometryHelper_t
Definition multidomain/facet/box/fvgridgeometry.hh:39
typename T::GeometryHelper SpecifiesGeometryHelper
Definition basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition method.hh:102
Definition adapt.hh:17
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
Definition common/pdesolver.hh:24
The default traits for the finite volume grid geometry of the box scheme with coupling occurring acro...
Definition multidomain/facet/box/fvgridgeometry.hh:55
BoxFacetCouplingFVElementGeometry< GridGeometry, enableCache > LocalView
Definition multidomain/facet/box/fvgridgeometry.hh:61
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition multidomain/facet/box/fvgridgeometry.hh:68
BoxFacetCouplingSubControlVolumeFace< GridView > SubControlVolumeFace
Definition multidomain/facet/box/fvgridgeometry.hh:58
EnrichedVertexDofMapper< GridView > VertexMapper
Definition multidomain/facet/box/fvgridgeometry.hh:66
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition multidomain/facet/box/fvgridgeometry.hh:64
BoxSubControlVolume< GridView > SubControlVolume
Definition multidomain/facet/box/fvgridgeometry.hh:57
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:27
unsigned int LocalIndex
Definition indextraits.hh:28
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...