3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
27#ifndef DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
28#define DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH
29
30#include <algorithm>
31#include <utility>
32
33#include <dune/grid/common/mcmgmapper.hh>
34#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
35
43
47
48namespace Dumux {
49
57template<class GridView>
59{
60 // use a specialized version of the box scvf
63
64 template<class GridGeometry, bool enableCache>
66
67 // per default we use an mcmg mapper for the elements
68 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
69 // the default vertex mapper is the enriched vertex dof mapper
71 // Mapper type for mapping edges
72 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
73};
74
82template<class Scalar,
83 class GridView,
84 bool enableGridGeometryCache = false,
87
95template<class Scalar, class GV, class Traits>
96class BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits>
97: public BaseGridGeometry<GV, Traits>
98{
101 using GridIndexType = typename IndexTraits<GV>::GridIndex;
102 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
103
104 using Element = typename GV::template Codim<0>::Entity;
105 using CoordScalar = typename GV::ctype;
106 static const int dim = GV::dimension;
107 static const int dimWorld = GV::dimensionworld;
108
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 [[deprecated("Use BoxFacetCouplingFVGridGeometry(gridView, facetGridView, codimOneGridAdapter) instead! Will be removed after release 3.5.")]]
134 : ParentType(gridView) {}
135
136 template<class FacetGridView, class CodimOneGridAdapter>
138 const FacetGridView& facetGridView,
139 const CodimOneGridAdapter& codimOneGridAdapter,
140 bool verbose = false)
141 : ParentType(gridView)
142 {
143 update_(facetGridView, codimOneGridAdapter, verbose);
144 }
145
147 const DofMapper& dofMapper() const
148 { return this->vertexMapper(); }
149
151 std::size_t numScv() const
152 { return numScv_; }
153
155 std::size_t numScvf() const
156 { return numScvf_; }
157
160 std::size_t numBoundaryScvf() const
161 { return numBoundaryScvf_; }
162
164 std::size_t numDofs() const
165 { return this->vertexMapper().size(); }
166
178 template<class FacetGridView, class CodimOneGridAdapter>
179 [[deprecated("Use update(gridView, facetGridView, codimOneGridAdapter) instead! Will be removed after release 3.5.")]]
180 void update(const FacetGridView& facetGridView,
181 const CodimOneGridAdapter& codimOneGridAdapter,
182 bool verbose = false)
183 {
184 // first update the parent (mappers etc)
185 ParentType::update();
186 update_(facetGridView, codimOneGridAdapter, verbose);
187 }
188
202 template<class FacetGridView, class CodimOneGridAdapter>
203 void update(const GridView& gridView,
204 const FacetGridView& facetGridView,
205 const CodimOneGridAdapter& codimOneGridAdapter,
206 bool verbose = false)
207 {
208 ParentType::update(gridView);
209 update_(facetGridView, codimOneGridAdapter, verbose);
210 }
211
213 template<class FacetGridView, class CodimOneGridAdapter>
214 void update(GridView&& gridView,
215 const FacetGridView& facetGridView,
216 const CodimOneGridAdapter& codimOneGridAdapter,
217 bool verbose = false)
218 {
219 ParentType::update(std::move(gridView));
220 update_(facetGridView, codimOneGridAdapter, verbose);
221 }
222
224 const FeCache& feCache() const
225 { return feCache_; }
226
228 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
229 { return scvs_[eIdx]; }
230
232 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
233 { return scvfs_[eIdx]; }
234
236 bool dofOnBoundary(GridIndexType dofIdx) const
237 { return boundaryDofIndices_[dofIdx]; }
238
240 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
241 { return interiorBoundaryDofIndices_[dofIdx]; }
242
244 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
245 { return false; }
246
248 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
249 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
250
252 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
253 { return std::unordered_map<GridIndexType, GridIndexType>(); }
254
255private:
256
257 template<class FacetGridView, class CodimOneGridAdapter>
258 void update_(const FacetGridView& facetGridView,
259 const CodimOneGridAdapter& codimOneGridAdapter,
260 bool verbose = false)
261 {
262 // enrich the vertex mapper subject to the provided facet grid
263 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
264
265 // resize containers
266 const auto numDof = numDofs();
267 const auto numElements = this->gridView().size(0);
268 scvs_.clear();
269 scvfs_.clear();
270 scvs_.resize(numElements);
271 scvfs_.resize(numElements);
272 boundaryDofIndices_.assign(numDof, false);
273 interiorBoundaryDofIndices_.assign(numDof, false);
274
275 // Build the SCV and SCV faces
276 numScv_ = 0;
277 numScvf_ = 0;
278 numBoundaryScvf_ = 0;
279 for (const auto& element : elements(this->gridView()))
280 {
281 auto eIdx = this->elementMapper().index(element);
282
283 // keep track of number of scvs and scvfs
284 numScv_ += element.subEntities(dim);
285 numScvf_ += element.subEntities(dim-1);
286
287 // get the element geometry
288 auto elementGeometry = element.geometry();
289 const auto refElement = referenceElement(elementGeometry);
290
291 // instantiate the geometry helper
292 GeometryHelper geometryHelper(elementGeometry);
293
294 // construct the sub control volumes
295 scvs_[eIdx].clear();
296 scvs_[eIdx].reserve(elementGeometry.corners());
297 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
298 scvs_[eIdx].emplace_back(geometryHelper,
299 scvLocalIdx,
300 eIdx,
301 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
302
303 // construct the sub control volume faces
304 LocalIndexType scvfLocalIdx = 0;
305 scvfs_[eIdx].clear();
306 scvfs_[eIdx].reserve(element.subEntities(dim-1));
307 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
308 {
309 // find the global and local scv indices this scvf is belonging to
310 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
311 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
312
313 // create the sub-control volume face
314 scvfs_[eIdx].emplace_back(geometryHelper,
315 element,
316 elementGeometry,
317 scvfLocalIdx,
318 std::move(localScvIndices));
319 }
320
321 // construct the sub control volume faces on the domain/interior boundaries
322 // skip handled facets (necessary for e.g. Dune::FoamGrid)
323 std::vector<unsigned int> handledFacets;
324 for (const auto& intersection : intersections(this->gridView(), element))
325 {
326 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
327 continue;
328
329 handledFacets.push_back(intersection.indexInInside());
330
331 // determine if all corners live on the facet grid
332 const auto isGeometry = intersection.geometry();
333 const auto numFaceCorners = isGeometry.corners();
334 const auto idxInInside = intersection.indexInInside();
335 const auto boundary = intersection.boundary();
336
337 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
338 for (int i = 0; i < numFaceCorners; ++i)
339 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
340
341 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
342 for (int i = 0; i < numFaceCorners; ++i)
343 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
344
345 // if the vertices compose a facet element, the intersection is on facet grid
346 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
347
348 // make sure there are no periodic boundaries
349 if (boundary && intersection.neighbor())
350 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
351
352 // if it is not, but it is on the boundary -> boundary scvf
353 if (isOnFacet || boundary)
354 {
355 // keep track of number of faces
356 numScvf_ += numFaceCorners;
357 numBoundaryScvf_ += int(boundary)*numFaceCorners;
358
359 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
360 {
361 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
362 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
363
364 // create the sub-control volume face
365 scvfs_[eIdx].emplace_back(geometryHelper,
366 intersection,
367 isGeometry,
368 isScvfLocalIdx,
369 scvfLocalIdx,
370 std::move(localScvIndices),
371 boundary,
372 isOnFacet);
373
374 // Mark vertices to be on domain and/or interior boundary
375 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
376 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary;
377 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
378
379 // increment local counter
380 scvfLocalIdx++;
381 }
382 }
383 }
384 }
385 }
386
387 const FeCache feCache_;
388
389 std::vector<std::vector<SubControlVolume>> scvs_;
390 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
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
409template<class Scalar, class GV, class Traits>
410class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
411: public BaseGridGeometry<GV, Traits>
412{
415 using GridIndexType = typename IndexTraits<GV>::GridIndex;
416 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
417
418 static const int dim = GV::dimension;
419 static const int dimWorld = GV::dimensionworld;
420
421 using Element = typename GV::template Codim<0>::Entity;
422 using Intersection = typename GV::Intersection;
423 using CoordScalar = typename GV::ctype;
424
425public:
428 static constexpr DiscretizationMethod discMethod{};
429
431 using LocalView = typename Traits::template LocalView<ThisType, false>;
433 using SubControlVolume = typename Traits::SubControlVolume;
435 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
439 using DofMapper = typename Traits::VertexMapper;
441 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
443 using GridView = GV;
444
446 [[deprecated("Use BoxFacetCouplingFVGridGeometry(gridView, facetGridView, codimOneGridAdapter) instead! Will be removed after release 3.5.")]]
448 : ParentType(gridView)
449 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
450 {}
451
452 template<class FacetGridView, class CodimOneGridAdapter>
454 const FacetGridView& facetGridView,
455 const CodimOneGridAdapter& codimOneGridAdapter,
456 bool verbose = false)
457 : ParentType(gridView)
458 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
459 {
460 update_(facetGridView, codimOneGridAdapter, verbose);
461 }
462
465 const DofMapper& dofMapper() const
466 { return this->vertexMapper(); }
467
469 std::size_t numScv() const
470 { return numScv_; }
471
473 std::size_t numScvf() const
474 { return numScvf_; }
475
478 std::size_t numBoundaryScvf() const
479 { return numBoundaryScvf_; }
480
482 std::size_t numDofs() const
483 { return this->vertexMapper().size(); }
484
496 template<class FacetGridView, class CodimOneGridAdapter>
497 [[deprecated("Use update(gridView, facetGridView, codimOneGridAdapter) instead! Will be removed after release 3.5.")]]
498 void update(const FacetGridView& facetGridView,
499 const CodimOneGridAdapter& codimOneGridAdapter,
500 bool verbose = false)
501 {
502 // first update the parent (mappers etc)
503 ParentType::update();
504 updateFacetMapper_();
505 update_(facetGridView, codimOneGridAdapter, verbose);
506 }
507
521 template<class FacetGridView, class CodimOneGridAdapter>
522 void update(const GridView& gridView,
523 const FacetGridView& facetGridView,
524 const CodimOneGridAdapter& codimOneGridAdapter,
525 bool verbose = false)
526 {
527 ParentType::update(gridView);
528 updateFacetMapper_();
529 update_(facetGridView, codimOneGridAdapter, verbose);
530 }
531
533 template<class FacetGridView, class CodimOneGridAdapter>
534 void update(GridView&& gridView,
535 const FacetGridView& facetGridView,
536 const CodimOneGridAdapter& codimOneGridAdapter,
537 bool verbose = false)
538 {
539 ParentType::update(std::move(gridView));
540 updateFacetMapper_();
541 update_(facetGridView, codimOneGridAdapter, verbose);
542 }
543
545 const FeCache& feCache() const
546 { return feCache_; }
547
549 bool dofOnBoundary(unsigned int dofIdx) const
550 { return boundaryDofIndices_[dofIdx]; }
551
553 bool dofOnInteriorBoundary(unsigned int dofIdx) const
554 { return interiorBoundaryDofIndices_[dofIdx]; }
555
557 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
558 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
559
561 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
562 { return false; }
563
565 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
566 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
567
569 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
570 { return std::unordered_map<GridIndexType, GridIndexType>(); }
571
572private:
573
574 void updateFacetMapper_()
575 {
576 if constexpr (Deprecated::hasUpdateGridView<typename Traits::FacetMapper, GridView>())
577 facetMapper_.update(this->gridView());
578 else
579 Deprecated::update(facetMapper_);
580 }
581
582 template<class FacetGridView, class CodimOneGridAdapter>
583 void update_(const FacetGridView& facetGridView,
584 const CodimOneGridAdapter& codimOneGridAdapter,
585 bool verbose)
586 {
587 // enrich the vertex mapper subject to the provided facet grid
588 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
589
590 // save global data on the grid's scvs and scvfs
591 // TODO do we need those information?
592 numScv_ = 0;
593 numScvf_ = 0;
594 numBoundaryScvf_ = 0;
595
596 const auto numDof = numDofs();
597 boundaryDofIndices_.assign(numDof, false);
598 interiorBoundaryDofIndices_.assign(numDof, false);
599 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
600 for (const auto& element : elements(this->gridView()))
601 {
602 numScv_ += element.subEntities(dim);
603 numScvf_ += element.subEntities(dim-1);
604
605 const auto elementGeometry = element.geometry();
606 const auto refElement = referenceElement(elementGeometry);
607
608 // store the sub control volume face indices on the domain/interior boundary
609 // skip handled facets (necessary for e.g. Dune::FoamGrid)
610 std::vector<unsigned int> handledFacets;
611 for (const auto& intersection : intersections(this->gridView(), element))
612 {
613 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
614 continue;
615
616 handledFacets.push_back(intersection.indexInInside());
617
618 // determine if all corners live on the facet grid
619 const auto isGeometry = intersection.geometry();
620 const auto numFaceCorners = isGeometry.corners();
621 const auto idxInInside = intersection.indexInInside();
622 const auto boundary = intersection.boundary();
623
624 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
625 for (int i = 0; i < numFaceCorners; ++i)
626 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
627
628 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
629 for (int i = 0; i < numFaceCorners; ++i)
630 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
631
632 // if all vertices are living on the facet grid, this is an interiour boundary
633 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
634
635 // make sure there are no periodic boundaries
636 if (boundary && intersection.neighbor())
637 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
638
639 if (isOnFacet || boundary)
640 {
641 numScvf_ += numFaceCorners;
642 numBoundaryScvf_ += int(boundary)*numFaceCorners;
643
644 // Mark vertices to be on domain and/or interior boundary
645 for (int i = 0; i < numFaceCorners; ++i)
646 {
647 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
648 if (boundary) boundaryDofIndices_[ dofIndex ] = true;
649 if (isOnFacet)
650 {
651 interiorBoundaryDofIndices_[ dofIndex ] = true;
652 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
653 }
654 }
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 domain/interior boundaries
669 std::vector<bool> boundaryDofIndices_;
670 std::vector<bool> interiorBoundaryDofIndices_;
671
672 // facet mapper and markers which facets lie on interior boundaries
673 typename Traits::FacetMapper facetMapper_;
674 std::vector<bool> facetIsOnInteriorBoundary_;
675};
676
677} // end namespace Dumux
678
679#endif
Defines the index types used for grid and local indices.
Helpers for deprecation.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Helper classes to compute the integration elements.
Base class for grid geometries.
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: common/pdesolver.hh:36
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all finite volume grid geometries.
Definition: basegridgeometry.hh:51
GV GridView
export the grid view type
Definition: basegridgeometry.hh:66
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:36
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:89
Definition: method.hh:73
Base class for the element-local finite volume geometry for box models in the context of models consi...
Definition: multidomain/facet/box/fvelementgeometry.hh:48
The default traits for the finite volume grid geometry of the box scheme with coupling occuring acros...
Definition: multidomain/facet/box/fvgridgeometry.hh:59
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:72
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:68
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:86
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:98
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:240
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:248
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:203
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:127
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:123
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:228
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:151
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:244
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:236
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:125
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:160
void update(const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (do this again after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:180
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:147
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:119
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:155
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:214
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:224
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:164
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Definition: multidomain/facet/box/fvgridgeometry.hh:137
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:117
BoxFacetCouplingFVGridGeometry(const GridView &gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:133
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:252
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:232
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:121
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:412
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:433
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Definition: multidomain/facet/box/fvgridgeometry.hh:453
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:469
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:473
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:482
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:569
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:534
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:478
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:441
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:545
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:439
void update(const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
update all fvElementGeometries (do this again after grid adaption)
Definition: multidomain/facet/box/fvgridgeometry.hh:498
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:465
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:431
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:565
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:561
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:435
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:437
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:549
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:522
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:553
BoxFacetCouplingFVGridGeometry(const GridView gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:447
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:557
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:52
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:52
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:192
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...
Definition: vertexmapper.hh:142
the sub control volume for the box scheme
Base class for a sub control volume face of the box method in the context of of models considering co...