version 3.10.0
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:
107 using DiscretizationMethod = DiscretizationMethods::Box;
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>
129 BoxFacetCouplingFVGridGeometry(const GridView& gridView,
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
221private:
222
223 template<class FacetGridView, class CodimOneGridAdapter>
224 void update_(const FacetGridView& facetGridView,
225 const CodimOneGridAdapter& codimOneGridAdapter,
226 bool verbose = false)
227 {
228 // enrich the vertex mapper subject to the provided facet grid
229 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
230
231 // resize containers
232 const auto numDof = numDofs();
233 const auto numElements = this->gridView().size(0);
234 scvs_.clear();
235 scvfs_.clear();
236 scvs_.resize(numElements);
237 scvfs_.resize(numElements);
238 boundaryDofIndices_.assign(numDof, false);
239 interiorBoundaryDofIndices_.assign(numDof, false);
240
241 // Build the SCV and SCV faces
242 numScv_ = 0;
243 numScvf_ = 0;
244 numBoundaryScvf_ = 0;
245 for (const auto& element : elements(this->gridView()))
246 {
247 auto eIdx = this->elementMapper().index(element);
248
249 // keep track of number of scvs and scvfs
250 numScv_ += element.subEntities(dim);
251 numScvf_ += element.subEntities(dim-1);
252
253 // get the element geometry
254 auto elementGeometry = element.geometry();
255 const auto refElement = referenceElement(elementGeometry);
256
257 // instantiate the geometry helper
258 GeometryHelper geometryHelper(elementGeometry);
259
260 // construct the sub control volumes
261 scvs_[eIdx].clear();
262 scvs_[eIdx].reserve(elementGeometry.corners());
263 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
264 scvs_[eIdx].emplace_back(geometryHelper.getScvCorners(scvLocalIdx),
265 scvLocalIdx,
266 eIdx,
267 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
268
269 // construct the sub control volume faces
270 LocalIndexType scvfLocalIdx = 0;
271 scvfs_[eIdx].clear();
272 scvfs_[eIdx].reserve(element.subEntities(dim-1));
273 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
274 {
275 // find the global and local scv indices this scvf is belonging to
276 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
277 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
278
279 // create the sub-control volume face
280 scvfs_[eIdx].emplace_back(geometryHelper,
281 element,
282 elementGeometry,
283 scvfLocalIdx,
284 std::move(localScvIndices));
285 }
286
287 // construct the sub control volume faces on the domain/interior boundaries
288 // skip handled facets (necessary for e.g. Dune::FoamGrid)
289 std::vector<unsigned int> handledFacets;
290 for (const auto& intersection : intersections(this->gridView(), element))
291 {
292 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
293 continue;
294
295 handledFacets.push_back(intersection.indexInInside());
296
297 // determine if all corners live on the facet grid
298 const auto isGeometry = intersection.geometry();
299 const auto numFaceCorners = isGeometry.corners();
300 const auto idxInInside = intersection.indexInInside();
301 const auto boundary = intersection.boundary();
302
303 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
304 for (int i = 0; i < numFaceCorners; ++i)
305 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
306
307 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
308 for (int i = 0; i < numFaceCorners; ++i)
309 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
310
311 // if the vertices compose a facet element, the intersection is on facet grid
312 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
313
314 // make sure there are no periodic boundaries
315 if (boundary && intersection.neighbor())
316 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
317
318 if (isOnFacet || boundary)
319 {
320 // keep track of number of faces
321 numScvf_ += numFaceCorners;
322 numBoundaryScvf_ += int(boundary)*numFaceCorners;
323
324 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
325 {
326 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
327 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
328
329 // create the sub-control volume face
330 scvfs_[eIdx].emplace_back(geometryHelper,
331 intersection,
332 isGeometry,
333 isScvfLocalIdx,
334 scvfLocalIdx,
335 std::move(localScvIndices),
336 boundary,
337 isOnFacet);
338
339 // Mark vertices to be on domain and/or interior boundary
340 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
341 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
342 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
343
344 // increment local counter
345 scvfLocalIdx++;
346 }
347 }
348 }
349 }
350 }
351
352 const FeCache feCache_;
353
354 std::vector<std::vector<SubControlVolume>> scvs_;
355 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
356
357 // TODO do we need those?
358 std::size_t numScv_;
359 std::size_t numScvf_;
360 std::size_t numBoundaryScvf_;
361
362 // vertices on domain/interior boundaries
363 std::vector<bool> boundaryDofIndices_;
364 std::vector<bool> interiorBoundaryDofIndices_;
365};
366
374template<class Scalar, class GV, class Traits>
375class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
376: public BaseGridGeometry<GV, Traits>
377{
378 using ThisType = BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>;
379 using ParentType = BaseGridGeometry<GV, Traits>;
380 using GridIndexType = typename IndexTraits<GV>::GridIndex;
381 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
382
383 static const int dim = GV::dimension;
384 static const int dimWorld = GV::dimensionworld;
385
386 using Element = typename GV::template Codim<0>::Entity;
387 using Intersection = typename GV::Intersection;
388 using CoordScalar = typename GV::ctype;
389
390public:
392 using DiscretizationMethod = DiscretizationMethods::Box;
393 static constexpr DiscretizationMethod discMethod{};
394
396 using LocalView = typename Traits::template LocalView<ThisType, false>;
398 using SubControlVolume = typename Traits::SubControlVolume;
400 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
402 using Extrusion = Extrusion_t<Traits>;
404 using DofMapper = typename Traits::VertexMapper;
406 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
408 using GridView = GV;
411
413 template<class FacetGridView, class CodimOneGridAdapter>
414 BoxFacetCouplingFVGridGeometry(const GridView& gridView,
415 const FacetGridView& facetGridView,
416 const CodimOneGridAdapter& codimOneGridAdapter,
417 bool verbose = false)
418 : ParentType(gridView)
419 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
420 {
421 update_(facetGridView, codimOneGridAdapter, verbose);
422 }
423
426 const DofMapper& dofMapper() const
427 { return this->vertexMapper(); }
428
430 std::size_t numScv() const
431 { return numScv_; }
432
434 std::size_t numScvf() const
435 { return numScvf_; }
436
439 std::size_t numBoundaryScvf() const
440 { return numBoundaryScvf_; }
441
443 std::size_t numDofs() const
444 { return this->vertexMapper().size(); }
445
459 template<class FacetGridView, class CodimOneGridAdapter>
460 void update(const GridView& gridView,
461 const FacetGridView& facetGridView,
462 const CodimOneGridAdapter& codimOneGridAdapter,
463 bool verbose = false)
464 {
465 ParentType::update(gridView);
466 updateFacetMapper_();
467 update_(facetGridView, codimOneGridAdapter, verbose);
468 }
469
471 template<class FacetGridView, class CodimOneGridAdapter>
472 void update(GridView&& gridView,
473 const FacetGridView& facetGridView,
474 const CodimOneGridAdapter& codimOneGridAdapter,
475 bool verbose = false)
476 {
477 ParentType::update(std::move(gridView));
478 updateFacetMapper_();
479 update_(facetGridView, codimOneGridAdapter, verbose);
480 }
481
483 const FeCache& feCache() const
484 { return feCache_; }
485
487 bool dofOnBoundary(unsigned int dofIdx) const
488 { return boundaryDofIndices_[dofIdx]; }
489
491 bool dofOnInteriorBoundary(unsigned int dofIdx) const
492 { return interiorBoundaryDofIndices_[dofIdx]; }
493
495 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
496 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
497
499 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
500 { return false; }
501
503 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
504 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
505
506private:
507
508 void updateFacetMapper_()
509 {
510 facetMapper_.update(this->gridView());
511 }
512
513 template<class FacetGridView, class CodimOneGridAdapter>
514 void update_(const FacetGridView& facetGridView,
515 const CodimOneGridAdapter& codimOneGridAdapter,
516 bool verbose)
517 {
518 // enrich the vertex mapper subject to the provided facet grid
519 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
520
521 // save global data on the grid's scvs and scvfs
522 // TODO do we need those information?
523 numScv_ = 0;
524 numScvf_ = 0;
525 numBoundaryScvf_ = 0;
526
527 const auto numDof = numDofs();
528 boundaryDofIndices_.assign(numDof, false);
529 interiorBoundaryDofIndices_.assign(numDof, false);
530 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
531 for (const auto& element : elements(this->gridView()))
532 {
533 numScv_ += element.subEntities(dim);
534 numScvf_ += element.subEntities(dim-1);
535
536 const auto elementGeometry = element.geometry();
537 const auto refElement = referenceElement(elementGeometry);
538
539 // store the sub control volume face indices on the domain/interior boundary
540 // skip handled facets (necessary for e.g. Dune::FoamGrid)
541 std::vector<unsigned int> handledFacets;
542 for (const auto& intersection : intersections(this->gridView(), element))
543 {
544 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
545 continue;
546
547 handledFacets.push_back(intersection.indexInInside());
548
549 // determine if all corners live on the facet grid
550 const auto isGeometry = intersection.geometry();
551 const auto numFaceCorners = isGeometry.corners();
552 const auto idxInInside = intersection.indexInInside();
553 const auto boundary = intersection.boundary();
554
555 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
556 for (int i = 0; i < numFaceCorners; ++i)
557 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
558
559 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
560 for (int i = 0; i < numFaceCorners; ++i)
561 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
562
563 // if all vertices are living on the facet grid, this is an interiour boundary
564 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
565
566 // make sure there are no periodic boundaries
567 if (boundary && intersection.neighbor())
568 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
569
570 if (isOnFacet || boundary)
571 {
572 numScvf_ += numFaceCorners;
573 numBoundaryScvf_ += int(boundary)*numFaceCorners;
574
575 // Mark vertices to be on domain and/or interior boundary
576 for (int i = 0; i < numFaceCorners; ++i)
577 {
578 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
579 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet;
580 if (isOnFacet)
581 {
582 interiorBoundaryDofIndices_[ dofIndex ] = true;
583 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
584 }
585 }
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 domain/interior boundaries
600 std::vector<bool> boundaryDofIndices_;
601 std::vector<bool> interiorBoundaryDofIndices_;
602
603 // facet mapper and markers which facets lie on interior boundaries
604 typename Traits::FacetMapper facetMapper_;
605 std::vector<bool> facetIsOnInteriorBoundary_;
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
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 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
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...
Definition vertexmapper.hh:126
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition extrusion.hh:155
the sub control volume for the box scheme
Helper classes to compute the integration elements.
@ element
Definition fieldtype.hh:23
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...
Distance implementation details.
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:94
Definition adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
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 ...