3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
27#define DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
28
29#include <utility>
30#include <unordered_map>
31
32#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
33
43
44namespace Dumux {
45
52template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
54: public MapperTraits
55{
58
59 template<class GridGeometry, bool enableCache>
61};
62
69template<class Scalar,
70 class GridView,
71 bool enableGridGeometryCache = false,
74
81template<class Scalar, class GV, class Traits>
82class BoxFVGridGeometry<Scalar, GV, true, Traits>
83: public BaseGridGeometry<GV, Traits>
84{
87 using GridIndexType = typename IndexTraits<GV>::GridIndex;
88 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
89
90 using Element = typename GV::template Codim<0>::Entity;
91 using CoordScalar = typename GV::ctype;
92 static const int dim = GV::dimension;
93 static const int dimWorld = GV::dimensionworld;
94
95 using GeometryHelper = BoxGeometryHelper<GV, dim,
96 typename Traits::SubControlVolume,
97 typename Traits::SubControlVolumeFace>;
98
99public:
102 static constexpr DiscretizationMethod discMethod{};
103
105 using LocalView = typename Traits::template LocalView<ThisType, true>;
107 using SubControlVolume = typename Traits::SubControlVolume;
109 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
113 using DofMapper = typename Traits::VertexMapper;
115 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
117 using GridView = GV;
118
121 : ParentType(gridView)
122 {
123 update_();
124 }
125
128 const DofMapper& dofMapper() const
129 { return this->vertexMapper(); }
130
132 std::size_t numScv() const
133 { return numScv_; }
134
136 std::size_t numScvf() const
137 { return numScvf_; }
138
141 std::size_t numBoundaryScvf() const
142 { return numBoundaryScvf_; }
143
145 std::size_t numDofs() const
146 { return this->vertexMapper().size(); }
147
149 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
150 void update()
151 {
152 ParentType::update();
153 update_();
154 }
155
157 void update(const GridView& gridView)
158 {
159 ParentType::update(gridView);
160 update_();
161 }
162
164 void update(GridView&& gridView)
165 {
166 ParentType::update(std::move(gridView));
167 update_();
168 }
169
171 const FeCache& feCache() const
172 { return feCache_; }
173
175 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
176 { return scvs_[eIdx]; }
177
179 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
180 { return scvfs_[eIdx]; }
181
183 bool dofOnBoundary(GridIndexType dofIdx) const
184 { return boundaryDofIndices_[dofIdx]; }
185
187 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
188 { return periodicVertexMap_.count(dofIdx); }
189
191 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
192 { return periodicVertexMap_.at(dofIdx); }
193
195 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
196 { return periodicVertexMap_; }
197
199 bool hasBoundaryScvf(GridIndexType eIdx) const
200 { return hasBoundaryScvf_[eIdx]; }
201
202private:
203 void update_()
204 {
205 scvs_.clear();
206 scvfs_.clear();
207
208 auto numElements = this->gridView().size(0);
209 scvs_.resize(numElements);
210 scvfs_.resize(numElements);
211 hasBoundaryScvf_.resize(numElements, false);
212
213 boundaryDofIndices_.assign(numDofs(), false);
214
215 numScv_ = 0;
216 numScvf_ = 0;
217 numBoundaryScvf_ = 0;
218 // Build the SCV and SCV faces
219 for (const auto& element : elements(this->gridView()))
220 {
221 // fill the element map with seeds
222 auto eIdx = this->elementMapper().index(element);
223
224 // count
225 numScv_ += element.subEntities(dim);
226 numScvf_ += element.subEntities(dim-1);
227
228 // get the element geometry
229 auto elementGeometry = element.geometry();
230 const auto refElement = referenceElement(elementGeometry);
231
232 // instantiate the geometry helper
233 GeometryHelper geometryHelper(elementGeometry);
234
235 // construct the sub control volumes
236 scvs_[eIdx].resize(elementGeometry.corners());
237 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
238 {
239 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
240
241 scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
242 scvLocalIdx,
243 eIdx,
244 dofIdxGlobal);
245 }
246
247 // construct the sub control volume faces
248 LocalIndexType scvfLocalIdx = 0;
249 scvfs_[eIdx].resize(element.subEntities(dim-1));
250 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
251 {
252 // find the global and local scv indices this scvf is belonging to
253 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
254 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
255
256 scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
257 element,
258 elementGeometry,
259 scvfLocalIdx,
260 std::move(localScvIndices),
261 false);
262 }
263
264 // construct the sub control volume faces on the domain boundary
265 for (const auto& intersection : intersections(this->gridView(), element))
266 {
267 if (intersection.boundary() && !intersection.neighbor())
268 {
269 const auto isGeometry = intersection.geometry();
270 hasBoundaryScvf_[eIdx] = true;
271
272 // count
273 numScvf_ += isGeometry.corners();
274 numBoundaryScvf_ += isGeometry.corners();
275
276 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
277 {
278 // find the scvs this scvf is belonging to
279 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
280 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
281
282 scvfs_[eIdx].emplace_back(geometryHelper,
283 intersection,
284 isGeometry,
285 isScvfLocalIdx,
286 scvfLocalIdx,
287 std::move(localScvIndices),
288 true);
289
290 // increment local counter
291 scvfLocalIdx++;
292 }
293
294 // add all vertices on the intersection to the set of
295 // boundary vertices
296 const auto fIdx = intersection.indexInInside();
297 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
298 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
299 {
300 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
301 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
302 boundaryDofIndices_[vIdxGlobal] = true;
303 }
304 }
305
306 // inform the grid geometry if we have periodic boundaries
307 else if (intersection.boundary() && intersection.neighbor())
308 {
309 this->setPeriodic();
310
311 // find the mapped periodic vertex of all vertices on periodic boundaries
312 const auto fIdx = intersection.indexInInside();
313 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
314 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
315 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
316 {
317 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
318 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
319 const auto vPos = elementGeometry.corner(vIdx);
320
321 const auto& outside = intersection.outside();
322 const auto outsideGeometry = outside.geometry();
323 for (const auto& isOutside : intersections(this->gridView(), outside))
324 {
325 // only check periodic vertices of the periodic neighbor
326 if (isOutside.boundary() && isOutside.neighbor())
327 {
328 const auto fIdxOutside = isOutside.indexInInside();
329 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
330 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
331 {
332 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
333 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
334 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
335 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
336 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
337 }
338 }
339 }
340 }
341 }
342 }
343 }
344
345 // error check: periodic boundaries currently don't work for box in parallel
346 if (this->isPeriodic() && this->gridView().comm().size() > 1)
347 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
348 }
349
350 const FeCache feCache_;
351
352 std::vector<std::vector<SubControlVolume>> scvs_;
353 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
354 // TODO do we need those?
355 std::size_t numScv_;
356 std::size_t numScvf_;
357 std::size_t numBoundaryScvf_;
358
359 // vertices on the boudary
360 std::vector<bool> boundaryDofIndices_;
361 std::vector<bool> hasBoundaryScvf_;
362
363 // a map for periodic boundary vertices
364 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
365};
366
374template<class Scalar, class GV, class Traits>
375class BoxFVGridGeometry<Scalar, GV, false, Traits>
376: public BaseGridGeometry<GV, Traits>
377{
380 using GridIndexType = typename IndexTraits<GV>::GridIndex;
381
382 static const int dim = GV::dimension;
383 static const int dimWorld = GV::dimensionworld;
384
385 using Element = typename GV::template Codim<0>::Entity;
386 using CoordScalar = typename GV::ctype;
387
388public:
391 static constexpr DiscretizationMethod discMethod{};
392
394 using LocalView = typename Traits::template LocalView<ThisType, false>;
396 using SubControlVolume = typename Traits::SubControlVolume;
398 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
402 using DofMapper = typename Traits::VertexMapper;
404 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
406 using GridView = GV;
407
410 : ParentType(gridView)
411 {
412 update_();
413 }
414
417 const DofMapper& dofMapper() const
418 { return this->vertexMapper(); }
419
421 std::size_t numScv() const
422 { return numScv_; }
423
425 std::size_t numScvf() const
426 { return numScvf_; }
427
430 std::size_t numBoundaryScvf() const
431 { return numBoundaryScvf_; }
432
434 std::size_t numDofs() const
435 { return this->vertexMapper().size(); }
436
438 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
439 void update()
440 {
441 ParentType::update();
442 update_();
443 }
444
446 void update(const GridView& gridView)
447 {
448 ParentType::update(gridView);
449 update_();
450 }
451
453 void update(GridView&& gridView)
454 {
455 ParentType::update(std::move(gridView));
456 update_();
457 }
458
460 const FeCache& feCache() const
461 { return feCache_; }
462
464 bool dofOnBoundary(GridIndexType dofIdx) const
465 { return boundaryDofIndices_[dofIdx]; }
466
468 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
469 { return periodicVertexMap_.count(dofIdx); }
470
472 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
473 { return periodicVertexMap_.at(dofIdx); }
474
476 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
477 { return periodicVertexMap_; }
478
479private:
480
481 void update_()
482 {
483 boundaryDofIndices_.assign(numDofs(), false);
484
485 // save global data on the grid's scvs and scvfs
486 // TODO do we need those information?
487 numScv_ = 0;
488 numScvf_ = 0;
489 numBoundaryScvf_ = 0;
490 for (const auto& element : elements(this->gridView()))
491 {
492 numScv_ += element.subEntities(dim);
493 numScvf_ += element.subEntities(dim-1);
494
495 const auto elementGeometry = element.geometry();
496 const auto refElement = referenceElement(elementGeometry);
497
498 // store the sub control volume face indices on the domain boundary
499 for (const auto& intersection : intersections(this->gridView(), element))
500 {
501 if (intersection.boundary() && !intersection.neighbor())
502 {
503 const auto isGeometry = intersection.geometry();
504 numScvf_ += isGeometry.corners();
505 numBoundaryScvf_ += isGeometry.corners();
506
507 // add all vertices on the intersection to the set of
508 // boundary vertices
509 const auto fIdx = intersection.indexInInside();
510 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
511 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
512 {
513 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
514 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
515 boundaryDofIndices_[vIdxGlobal] = true;
516 }
517 }
518
519 // inform the grid geometry if we have periodic boundaries
520 else if (intersection.boundary() && intersection.neighbor())
521 {
522 this->setPeriodic();
523
524 // find the mapped periodic vertex of all vertices on periodic boundaries
525 const auto fIdx = intersection.indexInInside();
526 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
527 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
528 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
529 {
530 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
531 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
532 const auto vPos = elementGeometry.corner(vIdx);
533
534 const auto& outside = intersection.outside();
535 const auto outsideGeometry = outside.geometry();
536 for (const auto& isOutside : intersections(this->gridView(), outside))
537 {
538 // only check periodic vertices of the periodic neighbor
539 if (isOutside.boundary() && isOutside.neighbor())
540 {
541 const auto fIdxOutside = isOutside.indexInInside();
542 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
543 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
544 {
545 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
546 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
547 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
548 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
549 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
550 }
551 }
552 }
553 }
554 }
555 }
556 }
557
558 // error check: periodic boundaries currently don't work for box in parallel
559 if (this->isPeriodic() && this->gridView().comm().size() > 1)
560 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
561 }
562
563 const FeCache feCache_;
564
565 // Information on the global number of geometries
566 // TODO do we need those information?
567 std::size_t numScv_;
568 std::size_t numScvf_;
569 std::size_t numBoundaryScvf_;
570
571 // vertices on the boudary
572 std::vector<bool> boundaryDofIndices_;
573
574 // a map for periodic boundary vertices
575 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
576};
577
578} // end namespace Dumux
579
580#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
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
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
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:49
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition: discretization/box/fvgridgeometry.hh:55
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:73
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:84
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: discretization/box/fvgridgeometry.hh:175
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: discretization/box/fvgridgeometry.hh:179
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:195
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:157
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:191
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvgridgeometry.hh:199
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:145
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:132
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:150
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:171
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:141
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:115
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:109
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:128
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:111
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:105
BoxFVGridGeometry(const GridView gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:120
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:183
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:136
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:187
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:164
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:107
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:113
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:377
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:430
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:468
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:394
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:398
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:439
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:396
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:472
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:460
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:425
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:402
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:464
BoxFVGridGeometry(const GridView gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:409
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:417
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:453
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:400
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:434
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:476
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:421
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:446
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:404
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:89
Class for a sub control volume face in the box method, i.e a part of the boundary of a sub control vo...
Definition: discretization/box/subcontrolvolumeface.hh:93
Definition: method.hh:73
Base class for the local finite volume geometry for box models This builds up the sub control volumes...
the sub control volume for the box scheme
Base class for a sub control volume face.