3.6-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
42
46
47namespace Dumux {
48
56template<class GridView>
58{
59 // use a specialized version of the box scvf
62
63 template<class GridGeometry, bool enableCache>
65
66 // per default we use an mcmg mapper for the elements
67 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
68 // the default vertex mapper is the enriched vertex dof mapper
70 // Mapper type for mapping edges
71 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
72};
73
81template<class Scalar,
82 class GridView,
83 bool enableGridGeometryCache = false,
86
94template<class Scalar, class GV, class Traits>
95class BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits>
96: public BaseGridGeometry<GV, Traits>
97{
100 using GridIndexType = typename IndexTraits<GV>::GridIndex;
101 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
102
103 using Element = typename GV::template Codim<0>::Entity;
104 using CoordScalar = typename GV::ctype;
105 static const int dim = GV::dimension;
106 static const int dimWorld = GV::dimensionworld;
107
109
110public:
113 static constexpr DiscretizationMethod discMethod{};
114
116 using LocalView = typename Traits::template LocalView<ThisType, true>;
118 using SubControlVolume = typename Traits::SubControlVolume;
120 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
124 using DofMapper = typename Traits::VertexMapper;
126 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
128 using GridView = GV;
129
131 template<class FacetGridView, class CodimOneGridAdapter>
133 const FacetGridView& facetGridView,
134 const CodimOneGridAdapter& codimOneGridAdapter,
135 bool verbose = false)
136 : ParentType(gridView)
137 {
138 update_(facetGridView, codimOneGridAdapter, verbose);
139 }
140
142 const DofMapper& dofMapper() const
143 { return this->vertexMapper(); }
144
146 std::size_t numScv() const
147 { return numScv_; }
148
150 std::size_t numScvf() const
151 { return numScvf_; }
152
155 std::size_t numBoundaryScvf() const
156 { return numBoundaryScvf_; }
157
159 std::size_t numDofs() const
160 { return this->vertexMapper().size(); }
161
175 template<class FacetGridView, class CodimOneGridAdapter>
176 void update(const GridView& gridView,
177 const FacetGridView& facetGridView,
178 const CodimOneGridAdapter& codimOneGridAdapter,
179 bool verbose = false)
180 {
181 ParentType::update(gridView);
182 update_(facetGridView, codimOneGridAdapter, verbose);
183 }
184
186 template<class FacetGridView, class CodimOneGridAdapter>
187 void update(GridView&& gridView,
188 const FacetGridView& facetGridView,
189 const CodimOneGridAdapter& codimOneGridAdapter,
190 bool verbose = false)
191 {
192 ParentType::update(std::move(gridView));
193 update_(facetGridView, codimOneGridAdapter, verbose);
194 }
195
197 const FeCache& feCache() const
198 { return feCache_; }
199
201 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
202 { return scvs_[eIdx]; }
203
205 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
206 { return scvfs_[eIdx]; }
207
209 bool dofOnBoundary(GridIndexType dofIdx) const
210 { return boundaryDofIndices_[dofIdx]; }
211
213 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
214 { return interiorBoundaryDofIndices_[dofIdx]; }
215
217 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
218 { return false; }
219
221 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
222 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
223
225 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
226 { return std::unordered_map<GridIndexType, GridIndexType>(); }
227
228private:
229
230 template<class FacetGridView, class CodimOneGridAdapter>
231 void update_(const FacetGridView& facetGridView,
232 const CodimOneGridAdapter& codimOneGridAdapter,
233 bool verbose = false)
234 {
235 // enrich the vertex mapper subject to the provided facet grid
236 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
237
238 // resize containers
239 const auto numDof = numDofs();
240 const auto numElements = this->gridView().size(0);
241 scvs_.clear();
242 scvfs_.clear();
243 scvs_.resize(numElements);
244 scvfs_.resize(numElements);
245 boundaryDofIndices_.assign(numDof, false);
246 interiorBoundaryDofIndices_.assign(numDof, false);
247
248 // Build the SCV and SCV faces
249 numScv_ = 0;
250 numScvf_ = 0;
251 numBoundaryScvf_ = 0;
252 for (const auto& element : elements(this->gridView()))
253 {
254 auto eIdx = this->elementMapper().index(element);
255
256 // keep track of number of scvs and scvfs
257 numScv_ += element.subEntities(dim);
258 numScvf_ += element.subEntities(dim-1);
259
260 // get the element geometry
261 auto elementGeometry = element.geometry();
262 const auto refElement = referenceElement(elementGeometry);
263
264 // instantiate the geometry helper
265 GeometryHelper geometryHelper(elementGeometry);
266
267 // construct the sub control volumes
268 scvs_[eIdx].clear();
269 scvs_[eIdx].reserve(elementGeometry.corners());
270 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
271 scvs_[eIdx].emplace_back(geometryHelper,
272 scvLocalIdx,
273 eIdx,
274 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
275
276 // construct the sub control volume faces
277 LocalIndexType scvfLocalIdx = 0;
278 scvfs_[eIdx].clear();
279 scvfs_[eIdx].reserve(element.subEntities(dim-1));
280 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
281 {
282 // find the global and local scv indices this scvf is belonging to
283 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
284 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
285
286 // create the sub-control volume face
287 scvfs_[eIdx].emplace_back(geometryHelper,
288 element,
289 elementGeometry,
290 scvfLocalIdx,
291 std::move(localScvIndices));
292 }
293
294 // construct the sub control volume faces on the domain/interior boundaries
295 // skip handled facets (necessary for e.g. Dune::FoamGrid)
296 std::vector<unsigned int> handledFacets;
297 for (const auto& intersection : intersections(this->gridView(), element))
298 {
299 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
300 continue;
301
302 handledFacets.push_back(intersection.indexInInside());
303
304 // determine if all corners live on the facet grid
305 const auto isGeometry = intersection.geometry();
306 const auto numFaceCorners = isGeometry.corners();
307 const auto idxInInside = intersection.indexInInside();
308 const auto boundary = intersection.boundary();
309
310 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
311 for (int i = 0; i < numFaceCorners; ++i)
312 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
313
314 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
315 for (int i = 0; i < numFaceCorners; ++i)
316 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
317
318 // if the vertices compose a facet element, the intersection is on facet grid
319 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
320
321 // make sure there are no periodic boundaries
322 if (boundary && intersection.neighbor())
323 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
324
325 // if it is not, but it is on the boundary -> boundary scvf
326 if (isOnFacet || boundary)
327 {
328 // keep track of number of faces
329 numScvf_ += numFaceCorners;
330 numBoundaryScvf_ += int(boundary)*numFaceCorners;
331
332 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
333 {
334 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
335 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
336
337 // create the sub-control volume face
338 scvfs_[eIdx].emplace_back(geometryHelper,
339 intersection,
340 isGeometry,
341 isScvfLocalIdx,
342 scvfLocalIdx,
343 std::move(localScvIndices),
344 boundary,
345 isOnFacet);
346
347 // Mark vertices to be on domain and/or interior boundary
348 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
349 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary;
350 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
351
352 // increment local counter
353 scvfLocalIdx++;
354 }
355 }
356 }
357 }
358 }
359
360 const FeCache feCache_;
361
362 std::vector<std::vector<SubControlVolume>> scvs_;
363 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
364
365 // TODO do we need those?
366 std::size_t numScv_;
367 std::size_t numScvf_;
368 std::size_t numBoundaryScvf_;
369
370 // vertices on domain/interior boundaries
371 std::vector<bool> boundaryDofIndices_;
372 std::vector<bool> interiorBoundaryDofIndices_;
373};
374
382template<class Scalar, class GV, class Traits>
383class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
384: public BaseGridGeometry<GV, Traits>
385{
388 using GridIndexType = typename IndexTraits<GV>::GridIndex;
389 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
390
391 static const int dim = GV::dimension;
392 static const int dimWorld = GV::dimensionworld;
393
394 using Element = typename GV::template Codim<0>::Entity;
395 using Intersection = typename GV::Intersection;
396 using CoordScalar = typename GV::ctype;
397
398public:
401 static constexpr DiscretizationMethod discMethod{};
402
404 using LocalView = typename Traits::template LocalView<ThisType, false>;
406 using SubControlVolume = typename Traits::SubControlVolume;
408 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
412 using DofMapper = typename Traits::VertexMapper;
414 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
416 using GridView = GV;
417
419 template<class FacetGridView, class CodimOneGridAdapter>
421 const FacetGridView& facetGridView,
422 const CodimOneGridAdapter& codimOneGridAdapter,
423 bool verbose = false)
424 : ParentType(gridView)
425 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
426 {
427 update_(facetGridView, codimOneGridAdapter, verbose);
428 }
429
432 const DofMapper& dofMapper() const
433 { return this->vertexMapper(); }
434
436 std::size_t numScv() const
437 { return numScv_; }
438
440 std::size_t numScvf() const
441 { return numScvf_; }
442
445 std::size_t numBoundaryScvf() const
446 { return numBoundaryScvf_; }
447
449 std::size_t numDofs() const
450 { return this->vertexMapper().size(); }
451
465 template<class FacetGridView, class CodimOneGridAdapter>
466 void update(const GridView& gridView,
467 const FacetGridView& facetGridView,
468 const CodimOneGridAdapter& codimOneGridAdapter,
469 bool verbose = false)
470 {
471 ParentType::update(gridView);
472 updateFacetMapper_();
473 update_(facetGridView, codimOneGridAdapter, verbose);
474 }
475
477 template<class FacetGridView, class CodimOneGridAdapter>
478 void update(GridView&& gridView,
479 const FacetGridView& facetGridView,
480 const CodimOneGridAdapter& codimOneGridAdapter,
481 bool verbose = false)
482 {
483 ParentType::update(std::move(gridView));
484 updateFacetMapper_();
485 update_(facetGridView, codimOneGridAdapter, verbose);
486 }
487
489 const FeCache& feCache() const
490 { return feCache_; }
491
493 bool dofOnBoundary(unsigned int dofIdx) const
494 { return boundaryDofIndices_[dofIdx]; }
495
497 bool dofOnInteriorBoundary(unsigned int dofIdx) const
498 { return interiorBoundaryDofIndices_[dofIdx]; }
499
501 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
502 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
503
505 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
506 { return false; }
507
509 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
510 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
511
513 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
514 { return std::unordered_map<GridIndexType, GridIndexType>(); }
515
516private:
517
518 void updateFacetMapper_()
519 {
520 facetMapper_.update(this->gridView());
521 }
522
523 template<class FacetGridView, class CodimOneGridAdapter>
524 void update_(const FacetGridView& facetGridView,
525 const CodimOneGridAdapter& codimOneGridAdapter,
526 bool verbose)
527 {
528 // enrich the vertex mapper subject to the provided facet grid
529 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
530
531 // save global data on the grid's scvs and scvfs
532 // TODO do we need those information?
533 numScv_ = 0;
534 numScvf_ = 0;
535 numBoundaryScvf_ = 0;
536
537 const auto numDof = numDofs();
538 boundaryDofIndices_.assign(numDof, false);
539 interiorBoundaryDofIndices_.assign(numDof, false);
540 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
541 for (const auto& element : elements(this->gridView()))
542 {
543 numScv_ += element.subEntities(dim);
544 numScvf_ += element.subEntities(dim-1);
545
546 const auto elementGeometry = element.geometry();
547 const auto refElement = referenceElement(elementGeometry);
548
549 // store the sub control volume face indices on the domain/interior boundary
550 // skip handled facets (necessary for e.g. Dune::FoamGrid)
551 std::vector<unsigned int> handledFacets;
552 for (const auto& intersection : intersections(this->gridView(), element))
553 {
554 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
555 continue;
556
557 handledFacets.push_back(intersection.indexInInside());
558
559 // determine if all corners live on the facet grid
560 const auto isGeometry = intersection.geometry();
561 const auto numFaceCorners = isGeometry.corners();
562 const auto idxInInside = intersection.indexInInside();
563 const auto boundary = intersection.boundary();
564
565 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
566 for (int i = 0; i < numFaceCorners; ++i)
567 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
568
569 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
570 for (int i = 0; i < numFaceCorners; ++i)
571 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
572
573 // if all vertices are living on the facet grid, this is an interiour boundary
574 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
575
576 // make sure there are no periodic boundaries
577 if (boundary && intersection.neighbor())
578 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
579
580 if (isOnFacet || boundary)
581 {
582 numScvf_ += numFaceCorners;
583 numBoundaryScvf_ += int(boundary)*numFaceCorners;
584
585 // Mark vertices to be on domain and/or interior boundary
586 for (int i = 0; i < numFaceCorners; ++i)
587 {
588 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
589 if (boundary) boundaryDofIndices_[ dofIndex ] = true;
590 if (isOnFacet)
591 {
592 interiorBoundaryDofIndices_[ dofIndex ] = true;
593 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
594 }
595 }
596 }
597 }
598 }
599 }
600
601 const FeCache feCache_;
602
603 // Information on the global number of geometries
604 // TODO do we need those information?
605 std::size_t numScv_;
606 std::size_t numScvf_;
607 std::size_t numBoundaryScvf_;
608
609 // vertices on domain/interior boundaries
610 std::vector<bool> boundaryDofIndices_;
611 std::vector<bool> interiorBoundaryDofIndices_;
612
613 // facet mapper and markers which facets lie on interior boundaries
614 typename Traits::FacetMapper facetMapper_;
615 std::vector<bool> facetIsOnInteriorBoundary_;
616};
617
618} // end namespace Dumux
619
620#endif
Defines the index types used for grid and local indices.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for grid geometries.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
Definition: deprecated.hh:149
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all grid geometries.
Definition: basegridgeometry.hh:61
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:69
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:72
Definition: method.hh:58
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 occurring acro...
Definition: multidomain/facet/box/fvgridgeometry.hh:58
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:71
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:67
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:85
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:97
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:213
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:221
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:176
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:126
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:122
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:201
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:146
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:217
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:209
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:124
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:155
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:142
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:118
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:150
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:187
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:197
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:159
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:132
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:116
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:225
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:205
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:120
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:385
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:406
BoxFacetCouplingFVGridGeometry(const GridView &gridView, const FacetGridView &facetGridView, const CodimOneGridAdapter &codimOneGridAdapter, bool verbose=false)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:420
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:436
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:440
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:449
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:513
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:478
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:445
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:414
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:489
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:412
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:432
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:404
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:509
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:505
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:408
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:410
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:493
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:466
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:497
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:501
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:53
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...