version 3.8
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-FileCopyrightInfo: 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{
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:
108 static constexpr DiscretizationMethod discMethod{};
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;
126
128 template<class FacetGridView, class CodimOneGridAdapter>
130 const FacetGridView& facetGridView,
131 const CodimOneGridAdapter& codimOneGridAdapter,
132 bool verbose = false)
133 : ParentType(gridView)
134 {
135 update_(facetGridView, codimOneGridAdapter, verbose);
136 }
137
139 const DofMapper& dofMapper() const
140 { return this->vertexMapper(); }
141
143 std::size_t numScv() const
144 { return numScv_; }
145
147 std::size_t numScvf() const
148 { return numScvf_; }
149
152 std::size_t numBoundaryScvf() const
153 { return numBoundaryScvf_; }
154
156 std::size_t numDofs() const
157 { return this->vertexMapper().size(); }
158
172 template<class FacetGridView, class CodimOneGridAdapter>
173 void update(const GridView& gridView,
174 const FacetGridView& facetGridView,
175 const CodimOneGridAdapter& codimOneGridAdapter,
176 bool verbose = false)
177 {
178 ParentType::update(gridView);
179 update_(facetGridView, codimOneGridAdapter, verbose);
180 }
181
183 template<class FacetGridView, class CodimOneGridAdapter>
184 void update(GridView&& gridView,
185 const FacetGridView& facetGridView,
186 const CodimOneGridAdapter& codimOneGridAdapter,
187 bool verbose = false)
188 {
189 ParentType::update(std::move(gridView));
190 update_(facetGridView, codimOneGridAdapter, verbose);
191 }
192
194 const FeCache& feCache() const
195 { return feCache_; }
196
198 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
199 { return scvs_[eIdx]; }
200
202 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
203 { return scvfs_[eIdx]; }
204
206 bool dofOnBoundary(GridIndexType dofIdx) const
207 { return boundaryDofIndices_[dofIdx]; }
208
210 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
211 { return interiorBoundaryDofIndices_[dofIdx]; }
212
214 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
215 { return false; }
216
218 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
219 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
220
222 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
223 { return std::unordered_map<GridIndexType, GridIndexType>(); }
224
225private:
226
227 template<class FacetGridView, class CodimOneGridAdapter>
228 void update_(const FacetGridView& facetGridView,
229 const CodimOneGridAdapter& codimOneGridAdapter,
230 bool verbose = false)
231 {
232 // enrich the vertex mapper subject to the provided facet grid
233 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
234
235 // resize containers
236 const auto numDof = numDofs();
237 const auto numElements = this->gridView().size(0);
238 scvs_.clear();
239 scvfs_.clear();
240 scvs_.resize(numElements);
241 scvfs_.resize(numElements);
242 boundaryDofIndices_.assign(numDof, false);
243 interiorBoundaryDofIndices_.assign(numDof, false);
244
245 // Build the SCV and SCV faces
246 numScv_ = 0;
247 numScvf_ = 0;
248 numBoundaryScvf_ = 0;
249 for (const auto& element : elements(this->gridView()))
250 {
251 auto eIdx = this->elementMapper().index(element);
252
253 // keep track of number of scvs and scvfs
254 numScv_ += element.subEntities(dim);
255 numScvf_ += element.subEntities(dim-1);
256
257 // get the element geometry
258 auto elementGeometry = element.geometry();
259 const auto refElement = referenceElement(elementGeometry);
260
261 // instantiate the geometry helper
262 GeometryHelper geometryHelper(elementGeometry);
263
264 // construct the sub control volumes
265 scvs_[eIdx].clear();
266 scvs_[eIdx].reserve(elementGeometry.corners());
267 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
268 scvs_[eIdx].emplace_back(geometryHelper.getScvCorners(scvLocalIdx),
269 scvLocalIdx,
270 eIdx,
271 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
272
273 // construct the sub control volume faces
274 LocalIndexType scvfLocalIdx = 0;
275 scvfs_[eIdx].clear();
276 scvfs_[eIdx].reserve(element.subEntities(dim-1));
277 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
278 {
279 // find the global and local scv indices this scvf is belonging to
280 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
281 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
282
283 // create the sub-control volume face
284 scvfs_[eIdx].emplace_back(geometryHelper,
285 element,
286 elementGeometry,
287 scvfLocalIdx,
288 std::move(localScvIndices));
289 }
290
291 // construct the sub control volume faces on the domain/interior boundaries
292 // skip handled facets (necessary for e.g. Dune::FoamGrid)
293 std::vector<unsigned int> handledFacets;
294 for (const auto& intersection : intersections(this->gridView(), element))
295 {
296 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
297 continue;
298
299 handledFacets.push_back(intersection.indexInInside());
300
301 // determine if all corners live on the facet grid
302 const auto isGeometry = intersection.geometry();
303 const auto numFaceCorners = isGeometry.corners();
304 const auto idxInInside = intersection.indexInInside();
305 const auto boundary = intersection.boundary();
306
307 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
308 for (int i = 0; i < numFaceCorners; ++i)
309 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
310
311 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
312 for (int i = 0; i < numFaceCorners; ++i)
313 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
314
315 // if the vertices compose a facet element, the intersection is on facet grid
316 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
317
318 // make sure there are no periodic boundaries
319 if (boundary && intersection.neighbor())
320 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
321
322 // if it is not, but it is on the boundary -> boundary scvf
323 if (isOnFacet || boundary)
324 {
325 // keep track of number of faces
326 numScvf_ += numFaceCorners;
327 numBoundaryScvf_ += int(boundary)*numFaceCorners;
328
329 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
330 {
331 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
332 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
333
334 // create the sub-control volume face
335 scvfs_[eIdx].emplace_back(geometryHelper,
336 intersection,
337 isGeometry,
338 isScvfLocalIdx,
339 scvfLocalIdx,
340 std::move(localScvIndices),
341 boundary,
342 isOnFacet);
343
344 // Mark vertices to be on domain and/or interior boundary
345 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
346 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary;
347 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
348
349 // increment local counter
350 scvfLocalIdx++;
351 }
352 }
353 }
354 }
355 }
356
357 const FeCache feCache_;
358
359 std::vector<std::vector<SubControlVolume>> scvs_;
360 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
361
362 // TODO do we need those?
363 std::size_t numScv_;
364 std::size_t numScvf_;
365 std::size_t numBoundaryScvf_;
366
367 // vertices on domain/interior boundaries
368 std::vector<bool> boundaryDofIndices_;
369 std::vector<bool> interiorBoundaryDofIndices_;
370};
371
379template<class Scalar, class GV, class Traits>
380class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
381: public BaseGridGeometry<GV, Traits>
382{
385 using GridIndexType = typename IndexTraits<GV>::GridIndex;
386 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
387
388 static const int dim = GV::dimension;
389 static const int dimWorld = GV::dimensionworld;
390
391 using Element = typename GV::template Codim<0>::Entity;
392 using Intersection = typename GV::Intersection;
393 using CoordScalar = typename GV::ctype;
394
395public:
398 static constexpr DiscretizationMethod discMethod{};
399
401 using LocalView = typename Traits::template LocalView<ThisType, false>;
403 using SubControlVolume = typename Traits::SubControlVolume;
405 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
409 using DofMapper = typename Traits::VertexMapper;
411 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
413 using GridView = GV;
416
418 template<class FacetGridView, class CodimOneGridAdapter>
420 const FacetGridView& facetGridView,
421 const CodimOneGridAdapter& codimOneGridAdapter,
422 bool verbose = false)
423 : ParentType(gridView)
424 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
425 {
426 update_(facetGridView, codimOneGridAdapter, verbose);
427 }
428
431 const DofMapper& dofMapper() const
432 { return this->vertexMapper(); }
433
435 std::size_t numScv() const
436 { return numScv_; }
437
439 std::size_t numScvf() const
440 { return numScvf_; }
441
444 std::size_t numBoundaryScvf() const
445 { return numBoundaryScvf_; }
446
448 std::size_t numDofs() const
449 { return this->vertexMapper().size(); }
450
464 template<class FacetGridView, class CodimOneGridAdapter>
465 void update(const GridView& gridView,
466 const FacetGridView& facetGridView,
467 const CodimOneGridAdapter& codimOneGridAdapter,
468 bool verbose = false)
469 {
470 ParentType::update(gridView);
471 updateFacetMapper_();
472 update_(facetGridView, codimOneGridAdapter, verbose);
473 }
474
476 template<class FacetGridView, class CodimOneGridAdapter>
477 void update(GridView&& gridView,
478 const FacetGridView& facetGridView,
479 const CodimOneGridAdapter& codimOneGridAdapter,
480 bool verbose = false)
481 {
482 ParentType::update(std::move(gridView));
483 updateFacetMapper_();
484 update_(facetGridView, codimOneGridAdapter, verbose);
485 }
486
488 const FeCache& feCache() const
489 { return feCache_; }
490
492 bool dofOnBoundary(unsigned int dofIdx) const
493 { return boundaryDofIndices_[dofIdx]; }
494
496 bool dofOnInteriorBoundary(unsigned int dofIdx) const
497 { return interiorBoundaryDofIndices_[dofIdx]; }
498
500 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
501 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
502
504 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
505 { return false; }
506
508 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
509 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
510
512 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
513 { return std::unordered_map<GridIndexType, GridIndexType>(); }
514
515private:
516
517 void updateFacetMapper_()
518 {
519 facetMapper_.update(this->gridView());
520 }
521
522 template<class FacetGridView, class CodimOneGridAdapter>
523 void update_(const FacetGridView& facetGridView,
524 const CodimOneGridAdapter& codimOneGridAdapter,
525 bool verbose)
526 {
527 // enrich the vertex mapper subject to the provided facet grid
528 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
529
530 // save global data on the grid's scvs and scvfs
531 // TODO do we need those information?
532 numScv_ = 0;
533 numScvf_ = 0;
534 numBoundaryScvf_ = 0;
535
536 const auto numDof = numDofs();
537 boundaryDofIndices_.assign(numDof, false);
538 interiorBoundaryDofIndices_.assign(numDof, false);
539 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
540 for (const auto& element : elements(this->gridView()))
541 {
542 numScv_ += element.subEntities(dim);
543 numScvf_ += element.subEntities(dim-1);
544
545 const auto elementGeometry = element.geometry();
546 const auto refElement = referenceElement(elementGeometry);
547
548 // store the sub control volume face indices on the domain/interior boundary
549 // skip handled facets (necessary for e.g. Dune::FoamGrid)
550 std::vector<unsigned int> handledFacets;
551 for (const auto& intersection : intersections(this->gridView(), element))
552 {
553 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
554 continue;
555
556 handledFacets.push_back(intersection.indexInInside());
557
558 // determine if all corners live on the facet grid
559 const auto isGeometry = intersection.geometry();
560 const auto numFaceCorners = isGeometry.corners();
561 const auto idxInInside = intersection.indexInInside();
562 const auto boundary = intersection.boundary();
563
564 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
565 for (int i = 0; i < numFaceCorners; ++i)
566 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
567
568 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
569 for (int i = 0; i < numFaceCorners; ++i)
570 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
571
572 // if all vertices are living on the facet grid, this is an interiour boundary
573 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
574
575 // make sure there are no periodic boundaries
576 if (boundary && intersection.neighbor())
577 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
578
579 if (isOnFacet || boundary)
580 {
581 numScvf_ += numFaceCorners;
582 numBoundaryScvf_ += int(boundary)*numFaceCorners;
583
584 // Mark vertices to be on domain and/or interior boundary
585 for (int i = 0; i < numFaceCorners; ++i)
586 {
587 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
588 if (boundary) boundaryDofIndices_[ dofIndex ] = true;
589 if (isOnFacet)
590 {
591 interiorBoundaryDofIndices_[ dofIndex ] = true;
592 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
593 }
594 }
595 }
596 }
597 }
598 }
599
600 const FeCache feCache_;
601
602 // Information on the global number of geometries
603 // TODO do we need those information?
604 std::size_t numScv_;
605 std::size_t numScvf_;
606 std::size_t numBoundaryScvf_;
607
608 // vertices on domain/interior boundaries
609 std::vector<bool> boundaryDofIndices_;
610 std::vector<bool> interiorBoundaryDofIndices_;
611
612 // facet mapper and markers which facets lie on interior boundaries
613 typename Traits::FacetMapper facetMapper_;
614 std::vector<bool> facetIsOnInteriorBoundary_;
615};
616
617} // end namespace Dumux
618
619#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 element-local finite volume geometry for box models in the context of models consi...
Definition: multidomain/facet/box/fvelementgeometry.hh:36
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:382
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:403
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:419
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:435
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:439
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:448
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:512
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:477
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:444
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:411
Detail::BoxFacetCouplingGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: multidomain/facet/box/fvgridgeometry.hh:415
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:488
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:409
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:431
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:401
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:508
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:504
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:405
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:407
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:492
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:465
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:496
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:500
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:94
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:210
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:218
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:173
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
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:198
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:143
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:214
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:206
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:152
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:139
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:113
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:147
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:184
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:194
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:156
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:129
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
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:222
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:202
Detail::BoxFacetCouplingGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: multidomain/facet/box/fvgridgeometry.hh:125
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:41
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:60
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:130
the sub control volume for the box scheme
Helper classes to compute the integration elements.
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Base class for a sub control volume face of the box method in the context of of models considering co...
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:43
typename T::GeometryHelper SpecifiesGeometryHelper
Definition: basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
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
The default traits for the finite volume grid geometry of the box scheme with coupling occurring acro...
Definition: multidomain/facet/box/fvgridgeometry.hh:55
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:68
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:64
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26