3.1-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
32#include <dune/grid/common/mcmgmapper.hh>
33#include <dune/geometry/referenceelements.hh>
34#include <dune/localfunctions/lagrange/pqkfactory.hh>
35
41
45
46namespace Dumux {
47
55template<class GridView>
57{
58 // use a specialized version of the box scvf
61
62 template<class GridGeometry, bool enableCache>
64
65 // per default we use an mcmg mapper for the elements
66 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
67 // the default vertex mapper is the enriched vertex dof mapper
69 // Mapper type for mapping edges
70 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
71};
72
80template<class Scalar,
81 class GridView,
82 bool enableGridGeometryCache = false,
85
93template<class Scalar, class GV, class Traits>
94class BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits>
95: public BaseGridGeometry<GV, Traits>
96{
99 using GridIndexType = typename IndexTraits<GV>::GridIndex;
100 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
101
102 using Element = typename GV::template Codim<0>::Entity;
103 using CoordScalar = typename GV::ctype;
104 static const int dim = GV::dimension;
105 static const int dimWorld = GV::dimensionworld;
106
107 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
109
110public:
112 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
113
115 using LocalView = typename Traits::template LocalView<ThisType, true>;
117 using SubControlVolume = typename Traits::SubControlVolume;
119 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
121 using DofMapper = typename Traits::VertexMapper;
123 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
125 using GridView = GV;
126
129 : ParentType(gridView) {}
130
132 const DofMapper& dofMapper() const
133 { return this->vertexMapper(); }
134
136 std::size_t numScv() const
137 { return numScv_; }
138
140 std::size_t numScvf() const
141 { return numScvf_; }
142
145 std::size_t numBoundaryScvf() const
146 { return numBoundaryScvf_; }
147
149 std::size_t numDofs() const
150 { return this->vertexMapper().size(); }
151
163 template<class FacetGridView, class CodimOneGridAdapter>
164 void update(const FacetGridView& facetGridView,
165 const CodimOneGridAdapter& codimOneGridAdapter,
166 bool verbose = false)
167 {
168 // first update the parent (mappers etc)
169 ParentType::update();
170
171 // enrich the vertex mapper subject to the provided facet grid
172 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
173
174 // resize containers
175 const auto numDof = numDofs();
176 const auto numElements = this->gridView().size(0);
177 scvs_.clear();
178 scvfs_.clear();
179 scvs_.resize(numElements);
180 scvfs_.resize(numElements);
181 boundaryDofIndices_.assign(numDof, false);
182 interiorBoundaryDofIndices_.assign(numDof, false);
183
184 // Build the SCV and SCV faces
185 numScv_ = 0;
186 numScvf_ = 0;
187 numBoundaryScvf_ = 0;
188 for (const auto& element : elements(this->gridView()))
189 {
190 auto eIdx = this->elementMapper().index(element);
191
192 // keep track of number of scvs and scvfs
193 numScv_ += element.subEntities(dim);
194 numScvf_ += element.subEntities(dim-1);
195
196 // get the element geometry
197 auto elementGeometry = element.geometry();
198 const auto referenceElement = ReferenceElements::general(elementGeometry.type());
199
200 // instantiate the geometry helper
201 GeometryHelper geometryHelper(elementGeometry);
202
203 // construct the sub control volumes
204 scvs_[eIdx].clear();
205 scvs_[eIdx].reserve(elementGeometry.corners());
206 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
207 scvs_[eIdx].emplace_back(geometryHelper,
208 scvLocalIdx,
209 eIdx,
210 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
211
212 // construct the sub control volume faces
213 LocalIndexType scvfLocalIdx = 0;
214 scvfs_[eIdx].clear();
215 scvfs_[eIdx].reserve(element.subEntities(dim-1));
216 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
217 {
218 // find the global and local scv indices this scvf is belonging to
219 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
220 static_cast<LocalIndexType>(referenceElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
221
222 // create the sub-control volume face
223 scvfs_[eIdx].emplace_back(geometryHelper,
224 element,
225 elementGeometry,
226 scvfLocalIdx,
227 std::move(localScvIndices));
228 }
229
230 // construct the sub control volume faces on the domain/interior boundaries
231 // skip handled facets (necessary for e.g. Dune::FoamGrid)
232 std::vector<unsigned int> handledFacets;
233 for (const auto& intersection : intersections(this->gridView(), element))
234 {
235 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
236 continue;
237
238 handledFacets.push_back(intersection.indexInInside());
239
240 // determine if all corners live on the facet grid
241 const auto isGeometry = intersection.geometry();
242 const auto numFaceCorners = isGeometry.corners();
243 const auto idxInInside = intersection.indexInInside();
244 const auto boundary = intersection.boundary();
245
246 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
247 for (int i = 0; i < numFaceCorners; ++i)
248 vIndicesLocal[i] = static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, i, dim));
249
250 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
251 for (int i = 0; i < numFaceCorners; ++i)
252 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
253
254 // if the vertices compose a facet element, the intersection is on facet grid
255 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
256
257 // make sure there are no periodic boundaries
258 if (boundary && intersection.neighbor())
259 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
260
261 // if it is not, but it is on the boundary -> boundary scvf
262 if (isOnFacet || boundary)
263 {
264 // keep track of number of faces
265 numScvf_ += numFaceCorners;
266 numBoundaryScvf_ += int(boundary)*numFaceCorners;
267
268 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
269 {
270 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
271 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
272
273 // create the sub-control volume face
274 scvfs_[eIdx].emplace_back(geometryHelper,
275 intersection,
276 isGeometry,
277 isScvfLocalIdx,
278 scvfLocalIdx,
279 std::move(localScvIndices),
280 boundary,
281 isOnFacet);
282
283 // Mark vertices to be on domain and/or interior boundary
284 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
285 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary;
286 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
287
288 // increment local counter
289 scvfLocalIdx++;
290 }
291 }
292 }
293 }
294 }
295
297 const FeCache& feCache() const
298 { return feCache_; }
299
301 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
302 { return scvs_[eIdx]; }
303
305 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
306 { return scvfs_[eIdx]; }
307
309 bool dofOnBoundary(GridIndexType dofIdx) const
310 { return boundaryDofIndices_[dofIdx]; }
311
313 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
314 { return interiorBoundaryDofIndices_[dofIdx]; }
315
317 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
318 { return false; }
319
321 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
322 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
323
325 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
326 { return std::unordered_map<GridIndexType, GridIndexType>(); }
327
328private:
329 const FeCache feCache_;
330
331 std::vector<std::vector<SubControlVolume>> scvs_;
332 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
333
334 // TODO do we need those?
335 std::size_t numScv_;
336 std::size_t numScvf_;
337 std::size_t numBoundaryScvf_;
338
339 // vertices on domain/interior boundaries
340 std::vector<bool> boundaryDofIndices_;
341 std::vector<bool> interiorBoundaryDofIndices_;
342};
343
351template<class Scalar, class GV, class Traits>
352class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
353: public BaseGridGeometry<GV, Traits>
354{
357 using GridIndexType = typename IndexTraits<GV>::GridIndex;
358 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
359
360 static const int dim = GV::dimension;
361 static const int dimWorld = GV::dimensionworld;
362
363 using Element = typename GV::template Codim<0>::Entity;
364 using Intersection = typename GV::Intersection;
365 using CoordScalar = typename GV::ctype;
366
367 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
368
369public:
371 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
372
374 using LocalView = typename Traits::template LocalView<ThisType, false>;
376 using SubControlVolume = typename Traits::SubControlVolume;
378 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
380 using DofMapper = typename Traits::VertexMapper;
382 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
384 using GridView = GV;
385
388 : ParentType(gridView)
389 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
390 {}
391
394 const DofMapper& dofMapper() const
395 { return this->vertexMapper(); }
396
398 std::size_t numScv() const
399 { return numScv_; }
400
402 std::size_t numScvf() const
403 { return numScvf_; }
404
407 std::size_t numBoundaryScvf() const
408 { return numBoundaryScvf_; }
409
411 std::size_t numDofs() const
412 { return this->vertexMapper().size(); }
413
425 template<class FacetGridView, class CodimOneGridAdapter>
426 void update(const FacetGridView& facetGridView,
427 const CodimOneGridAdapter& codimOneGridAdapter,
428 bool verbose = false)
429 {
430 // first update the parent (mappers etc)
431 ParentType::update();
432
433 // enrich the vertex mapper subject to the provided facet grid
434 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
435
436 // save global data on the grid's scvs and scvfs
437 // TODO do we need those information?
438 numScv_ = 0;
439 numScvf_ = 0;
440 numBoundaryScvf_ = 0;
441
442 const auto numDof = numDofs();
443 boundaryDofIndices_.assign(numDof, false);
444 interiorBoundaryDofIndices_.assign(numDof, false);
445 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
446 for (const auto& element : elements(this->gridView()))
447 {
448 numScv_ += element.subEntities(dim);
449 numScvf_ += element.subEntities(dim-1);
450
451 const auto elementGeometry = element.geometry();
452 const auto referenceElement = ReferenceElements::general(elementGeometry.type());
453
454 // store the sub control volume face indices on the domain/interior boundary
455 // skip handled facets (necessary for e.g. Dune::FoamGrid)
456 std::vector<unsigned int> handledFacets;
457 for (const auto& intersection : intersections(this->gridView(), element))
458 {
459 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
460 continue;
461
462 handledFacets.push_back(intersection.indexInInside());
463
464 // determine if all corners live on the facet grid
465 const auto isGeometry = intersection.geometry();
466 const auto numFaceCorners = isGeometry.corners();
467 const auto idxInInside = intersection.indexInInside();
468 const auto boundary = intersection.boundary();
469
470 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
471 for (int i = 0; i < numFaceCorners; ++i)
472 vIndicesLocal[i] = static_cast<LocalIndexType>(referenceElement.subEntity(idxInInside, 1, i, dim));
473
474 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
475 for (int i = 0; i < numFaceCorners; ++i)
476 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
477
478 // if all vertices are living on the facet grid, this is an interiour boundary
479 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
480
481 // make sure there are no periodic boundaries
482 if (boundary && intersection.neighbor())
483 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
484
485 if (isOnFacet || boundary)
486 {
487 numScvf_ += numFaceCorners;
488 numBoundaryScvf_ += int(boundary)*numFaceCorners;
489
490 // Mark vertices to be on domain and/or interior boundary
491 for (int i = 0; i < numFaceCorners; ++i)
492 {
493 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
494 if (boundary) boundaryDofIndices_[ dofIndex ] = true;
495 if (isOnFacet)
496 {
497 interiorBoundaryDofIndices_[ dofIndex ] = true;
498 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
499 }
500 }
501 }
502 }
503 }
504 }
505
507 const FeCache& feCache() const
508 { return feCache_; }
509
511 bool dofOnBoundary(unsigned int dofIdx) const
512 { return boundaryDofIndices_[dofIdx]; }
513
515 bool dofOnInteriorBoundary(unsigned int dofIdx) const
516 { return interiorBoundaryDofIndices_[dofIdx]; }
517
519 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
520 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
521
523 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
524 { return false; }
525
527 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
528 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
529
531 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
532 { return std::unordered_map<GridIndexType, GridIndexType>(); }
533
534private:
535 const FeCache feCache_;
536
537 // Information on the global number of geometries
538 // TODO do we need those information?
539 std::size_t numScv_;
540 std::size_t numScvf_;
541 std::size_t numBoundaryScvf_;
542
543 // vertices on domain/interior boundaries
544 std::vector<bool> boundaryDofIndices_;
545 std::vector<bool> interiorBoundaryDofIndices_;
546
547 // facet mapper and markers which facets lie on interior boundaries
548 typename Traits::FacetMapper facetMapper_;
549 std::vector<bool> facetIsOnInteriorBoundary_;
550};
551
552} // end namespace Dumux
553
554#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.
The available discretization methods in Dumux.
copydoc Dumux::EnrichedVertexDofMapper
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
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:49
GV GridView
export the grid view type
Definition: basegridgeometry.hh:64
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:37
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:88
Base class for the element-local finite volume geometry for box models in the context of models consi...
Definition: multidomain/facet/box/fvelementgeometry.hh:48
The default traits for the finite volume grid geometry of the box scheme with coupling occuring acros...
Definition: multidomain/facet/box/fvgridgeometry.hh:57
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:70
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > ElementMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:66
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:84
Base class for the finite volume geometry vector for box schemes in the context of coupled models whe...
Definition: multidomain/facet/box/fvgridgeometry.hh:96
bool dofOnInteriorBoundary(GridIndexType dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:313
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:321
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:301
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:136
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:317
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:309
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:121
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:145
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:164
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:132
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:117
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:140
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:297
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:149
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:115
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:123
BoxFacetCouplingFVGridGeometry(const GridView &gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:128
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:325
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:305
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:119
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:354
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:376
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:398
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:402
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:411
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:531
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:407
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:507
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:380
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:426
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:394
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:374
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:527
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:523
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:382
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:378
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:511
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:515
BoxFacetCouplingFVGridGeometry(const GridView gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:387
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:519
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:193
A vertex mapper that allows for enrichment of nodes. Indication on where to enrich the nodes is done ...
Definition: vertexmapper.hh:146
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...