3.3.0
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/localfunctions/lagrange/pqkfactory.hh>
34
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
108
109public:
111 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
112
114 using LocalView = typename Traits::template LocalView<ThisType, true>;
116 using SubControlVolume = typename Traits::SubControlVolume;
118 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
122 using DofMapper = typename Traits::VertexMapper;
124 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
126 using GridView = GV;
127
130 : ParentType(gridView) {}
131
133 const DofMapper& dofMapper() const
134 { return this->vertexMapper(); }
135
137 std::size_t numScv() const
138 { return numScv_; }
139
141 std::size_t numScvf() const
142 { return numScvf_; }
143
146 std::size_t numBoundaryScvf() const
147 { return numBoundaryScvf_; }
148
150 std::size_t numDofs() const
151 { return this->vertexMapper().size(); }
152
164 template<class FacetGridView, class CodimOneGridAdapter>
165 void update(const FacetGridView& facetGridView,
166 const CodimOneGridAdapter& codimOneGridAdapter,
167 bool verbose = false)
168 {
169 // first update the parent (mappers etc)
170 ParentType::update();
171
172 // enrich the vertex mapper subject to the provided facet grid
173 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
174
175 // resize containers
176 const auto numDof = numDofs();
177 const auto numElements = this->gridView().size(0);
178 scvs_.clear();
179 scvfs_.clear();
180 scvs_.resize(numElements);
181 scvfs_.resize(numElements);
182 boundaryDofIndices_.assign(numDof, false);
183 interiorBoundaryDofIndices_.assign(numDof, false);
184
185 // Build the SCV and SCV faces
186 numScv_ = 0;
187 numScvf_ = 0;
188 numBoundaryScvf_ = 0;
189 for (const auto& element : elements(this->gridView()))
190 {
191 auto eIdx = this->elementMapper().index(element);
192
193 // keep track of number of scvs and scvfs
194 numScv_ += element.subEntities(dim);
195 numScvf_ += element.subEntities(dim-1);
196
197 // get the element geometry
198 auto elementGeometry = element.geometry();
199 const auto refElement = referenceElement(elementGeometry);
200
201 // instantiate the geometry helper
202 GeometryHelper geometryHelper(elementGeometry);
203
204 // construct the sub control volumes
205 scvs_[eIdx].clear();
206 scvs_[eIdx].reserve(elementGeometry.corners());
207 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
208 scvs_[eIdx].emplace_back(geometryHelper,
209 scvLocalIdx,
210 eIdx,
211 this->vertexMapper().subIndex(element, scvLocalIdx, dim));
212
213 // construct the sub control volume faces
214 LocalIndexType scvfLocalIdx = 0;
215 scvfs_[eIdx].clear();
216 scvfs_[eIdx].reserve(element.subEntities(dim-1));
217 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
218 {
219 // find the global and local scv indices this scvf is belonging to
220 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
221 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
222
223 // create the sub-control volume face
224 scvfs_[eIdx].emplace_back(geometryHelper,
225 element,
226 elementGeometry,
227 scvfLocalIdx,
228 std::move(localScvIndices));
229 }
230
231 // construct the sub control volume faces on the domain/interior boundaries
232 // skip handled facets (necessary for e.g. Dune::FoamGrid)
233 std::vector<unsigned int> handledFacets;
234 for (const auto& intersection : intersections(this->gridView(), element))
235 {
236 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
237 continue;
238
239 handledFacets.push_back(intersection.indexInInside());
240
241 // determine if all corners live on the facet grid
242 const auto isGeometry = intersection.geometry();
243 const auto numFaceCorners = isGeometry.corners();
244 const auto idxInInside = intersection.indexInInside();
245 const auto boundary = intersection.boundary();
246
247 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
248 for (int i = 0; i < numFaceCorners; ++i)
249 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
250
251 std::vector<LocalIndexType> gridVertexIndices(numFaceCorners);
252 for (int i = 0; i < numFaceCorners; ++i)
253 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
254
255 // if the vertices compose a facet element, the intersection is on facet grid
256 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
257
258 // make sure there are no periodic boundaries
259 if (boundary && intersection.neighbor())
260 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
261
262 // if it is not, but it is on the boundary -> boundary scvf
263 if (isOnFacet || boundary)
264 {
265 // keep track of number of faces
266 numScvf_ += numFaceCorners;
267 numBoundaryScvf_ += int(boundary)*numFaceCorners;
268
269 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx)
270 {
271 // find the inside scv this scvf is belonging to (localIdx = element local vertex index)
272 std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]};
273
274 // create the sub-control volume face
275 scvfs_[eIdx].emplace_back(geometryHelper,
276 intersection,
277 isGeometry,
278 isScvfLocalIdx,
279 scvfLocalIdx,
280 std::move(localScvIndices),
281 boundary,
282 isOnFacet);
283
284 // Mark vertices to be on domain and/or interior boundary
285 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim);
286 if (boundary) boundaryDofIndices_[ dofIndex ] = boundary;
287 if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet;
288
289 // increment local counter
290 scvfLocalIdx++;
291 }
292 }
293 }
294 }
295 }
296
298 const FeCache& feCache() const
299 { return feCache_; }
300
302 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
303 { return scvs_[eIdx]; }
304
306 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
307 { return scvfs_[eIdx]; }
308
310 bool dofOnBoundary(GridIndexType dofIdx) const
311 { return boundaryDofIndices_[dofIdx]; }
312
314 bool dofOnInteriorBoundary(GridIndexType dofIdx) const
315 { return interiorBoundaryDofIndices_[dofIdx]; }
316
318 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
319 { return false; }
320
322 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
323 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); }
324
326 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
327 { return std::unordered_map<GridIndexType, GridIndexType>(); }
328
329private:
330 const FeCache feCache_;
331
332 std::vector<std::vector<SubControlVolume>> scvs_;
333 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
334
335 // TODO do we need those?
336 std::size_t numScv_;
337 std::size_t numScvf_;
338 std::size_t numBoundaryScvf_;
339
340 // vertices on domain/interior boundaries
341 std::vector<bool> boundaryDofIndices_;
342 std::vector<bool> interiorBoundaryDofIndices_;
343};
344
352template<class Scalar, class GV, class Traits>
353class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>
354: public BaseGridGeometry<GV, Traits>
355{
358 using GridIndexType = typename IndexTraits<GV>::GridIndex;
359 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
360
361 static const int dim = GV::dimension;
362 static const int dimWorld = GV::dimensionworld;
363
364 using Element = typename GV::template Codim<0>::Entity;
365 using Intersection = typename GV::Intersection;
366 using CoordScalar = typename GV::ctype;
367
368public:
370 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
371
373 using LocalView = typename Traits::template LocalView<ThisType, false>;
375 using SubControlVolume = typename Traits::SubControlVolume;
377 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
381 using DofMapper = typename Traits::VertexMapper;
383 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
385 using GridView = GV;
386
389 : ParentType(gridView)
390 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
391 {}
392
395 const DofMapper& dofMapper() const
396 { return this->vertexMapper(); }
397
399 std::size_t numScv() const
400 { return numScv_; }
401
403 std::size_t numScvf() const
404 { return numScvf_; }
405
408 std::size_t numBoundaryScvf() const
409 { return numBoundaryScvf_; }
410
412 std::size_t numDofs() const
413 { return this->vertexMapper().size(); }
414
426 template<class FacetGridView, class CodimOneGridAdapter>
427 void update(const FacetGridView& facetGridView,
428 const CodimOneGridAdapter& codimOneGridAdapter,
429 bool verbose = false)
430 {
431 // first update the parent (mappers etc)
432 ParentType::update();
433
434 // enrich the vertex mapper subject to the provided facet grid
435 this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose);
436
437 // save global data on the grid's scvs and scvfs
438 // TODO do we need those information?
439 numScv_ = 0;
440 numScvf_ = 0;
441 numBoundaryScvf_ = 0;
442
443 const auto numDof = numDofs();
444 boundaryDofIndices_.assign(numDof, false);
445 interiorBoundaryDofIndices_.assign(numDof, false);
446 facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false);
447 for (const auto& element : elements(this->gridView()))
448 {
449 numScv_ += element.subEntities(dim);
450 numScvf_ += element.subEntities(dim-1);
451
452 const auto elementGeometry = element.geometry();
453 const auto refElement = referenceElement(elementGeometry);
454
455 // store the sub control volume face indices on the domain/interior boundary
456 // skip handled facets (necessary for e.g. Dune::FoamGrid)
457 std::vector<unsigned int> handledFacets;
458 for (const auto& intersection : intersections(this->gridView(), element))
459 {
460 if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside()))
461 continue;
462
463 handledFacets.push_back(intersection.indexInInside());
464
465 // determine if all corners live on the facet grid
466 const auto isGeometry = intersection.geometry();
467 const auto numFaceCorners = isGeometry.corners();
468 const auto idxInInside = intersection.indexInInside();
469 const auto boundary = intersection.boundary();
470
471 std::vector<LocalIndexType> vIndicesLocal(numFaceCorners);
472 for (int i = 0; i < numFaceCorners; ++i)
473 vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim));
474
475 std::vector<GridIndexType> gridVertexIndices(numFaceCorners);
476 for (int i = 0; i < numFaceCorners; ++i)
477 gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim);
478
479 // if all vertices are living on the facet grid, this is an interiour boundary
480 const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices);
481
482 // make sure there are no periodic boundaries
483 if (boundary && intersection.neighbor())
484 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme");
485
486 if (isOnFacet || boundary)
487 {
488 numScvf_ += numFaceCorners;
489 numBoundaryScvf_ += int(boundary)*numFaceCorners;
490
491 // Mark vertices to be on domain and/or interior boundary
492 for (int i = 0; i < numFaceCorners; ++i)
493 {
494 const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim);
495 if (boundary) boundaryDofIndices_[ dofIndex ] = true;
496 if (isOnFacet)
497 {
498 interiorBoundaryDofIndices_[ dofIndex ] = true;
499 facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true;
500 }
501 }
502 }
503 }
504 }
505 }
506
508 const FeCache& feCache() const
509 { return feCache_; }
510
512 bool dofOnBoundary(unsigned int dofIdx) const
513 { return boundaryDofIndices_[dofIdx]; }
514
516 bool dofOnInteriorBoundary(unsigned int dofIdx) const
517 { return interiorBoundaryDofIndices_[dofIdx]; }
518
520 bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const
521 { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; }
522
524 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
525 { return false; }
526
528 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
529 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); }
530
532 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
533 { return std::unordered_map<GridIndexType, GridIndexType>(); }
534
535private:
536 const FeCache feCache_;
537
538 // Information on the global number of geometries
539 // TODO do we need those information?
540 std::size_t numScv_;
541 std::size_t numScvf_;
542 std::size_t numBoundaryScvf_;
543
544 // vertices on domain/interior boundaries
545 std::vector<bool> boundaryDofIndices_;
546 std::vector<bool> interiorBoundaryDofIndices_;
547
548 // facet mapper and markers which facets lie on interior boundaries
549 typename Traits::FacetMapper facetMapper_;
550 std::vector<bool> facetIsOnInteriorBoundary_;
551};
552
553} // end namespace Dumux
554
555#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.
The available discretization methods in Dumux.
Base class for grid geometries.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:35
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:36
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:47
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:314
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:322
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:120
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:302
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:137
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:318
bool dofOnBoundary(GridIndexType dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:310
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:122
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:146
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:165
const DofMapper & dofMapper() const
the vertex mapper is the dofMapper
Definition: multidomain/facet/box/fvgridgeometry.hh:133
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:116
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:141
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:298
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:150
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:114
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:124
BoxFacetCouplingFVGridGeometry(const GridView &gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:129
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:326
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: multidomain/facet/box/fvgridgeometry.hh:306
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:118
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: multidomain/facet/box/fvgridgeometry.hh:355
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:375
std::size_t numScv() const
The total number of sub control volumes.
Definition: multidomain/facet/box/fvgridgeometry.hh:399
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: multidomain/facet/box/fvgridgeometry.hh:403
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: multidomain/facet/box/fvgridgeometry.hh:412
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: multidomain/facet/box/fvgridgeometry.hh:532
std::size_t numBoundaryScvf() const
Definition: multidomain/facet/box/fvgridgeometry.hh:408
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: multidomain/facet/box/fvgridgeometry.hh:508
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: multidomain/facet/box/fvgridgeometry.hh:381
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:427
const DofMapper & dofMapper() const
Definition: multidomain/facet/box/fvgridgeometry.hh:395
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:373
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:528
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
Periodic boundaries are not supported for the box facet coupling scheme.
Definition: multidomain/facet/box/fvgridgeometry.hh:524
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: multidomain/facet/box/fvgridgeometry.hh:383
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: multidomain/facet/box/fvgridgeometry.hh:377
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: multidomain/facet/box/fvgridgeometry.hh:379
bool dofOnBoundary(unsigned int dofIdx) const
If a d.o.f. is on the boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:512
bool dofOnInteriorBoundary(unsigned int dofIdx) const
If a d.o.f. is on an interior boundary.
Definition: multidomain/facet/box/fvgridgeometry.hh:516
BoxFacetCouplingFVGridGeometry(const GridView gridView)
Constructor.
Definition: multidomain/facet/box/fvgridgeometry.hh:388
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:520
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...