version 3.9
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
16
17#include <utility>
18#include <unordered_map>
19#include <array>
20#include <vector>
21
22#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
23
33
34namespace Dumux {
35
36namespace Detail {
37template<class GV, class T>
38using BoxGeometryHelper_t = Dune::Std::detected_or_t<
41 T
42>;
43} // end namespace Detail
44
51template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
53: public MapperTraits
54{
57
58 template<class GridGeometry, bool enableCache>
60};
61
68template<class Scalar,
69 class GridView,
70 bool enableGridGeometryCache = false,
73
80template<class Scalar, class GV, class Traits>
81class BoxFVGridGeometry<Scalar, GV, true, Traits>
82: public BaseGridGeometry<GV, Traits>
83{
86 using GridIndexType = typename IndexTraits<GV>::GridIndex;
87 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
88
89 using Element = typename GV::template Codim<0>::Entity;
90 using CoordScalar = typename GV::ctype;
91 static const int dim = GV::dimension;
92 static const int dimWorld = GV::dimensionworld;
93
94public:
97 static constexpr DiscretizationMethod discMethod{};
98
102 using LocalView = typename Traits::template LocalView<ThisType, true>;
104 using SubControlVolume = typename Traits::SubControlVolume;
106 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
110 using DofMapper = typename Traits::VertexMapper;
112 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
114 using GridView = GV;
115
117 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
118 : ParentType(std::move(gg))
119 , cache_(*this)
120 {
121 update_();
122 }
123
126 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
127 {}
128
131 const DofMapper& dofMapper() const
132 { return this->vertexMapper(); }
133
135 std::size_t numScv() const
136 { return numScv_; }
137
139 std::size_t numScvf() const
140 { return numScvf_; }
141
144 std::size_t numBoundaryScvf() const
145 { return numBoundaryScvf_; }
146
148 std::size_t numDofs() const
149 { return this->vertexMapper().size(); }
150
151
153 void update(const GridView& gridView)
154 {
155 ParentType::update(gridView);
156 update_();
157 }
158
160 void update(GridView&& gridView)
161 {
162 ParentType::update(std::move(gridView));
163 update_();
164 }
165
167 const FeCache& feCache() const
168 { return feCache_; }
169
171 bool dofOnBoundary(GridIndexType dofIdx) const
172 { return boundaryDofIndices_[dofIdx]; }
173
175 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
176 { return periodicVertexMap_.count(dofIdx); }
177
179 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
180 { return periodicVertexMap_.at(dofIdx); }
181
183 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
184 { return periodicVertexMap_; }
185
187 [[deprecated("Will be removed after release 3.9. Use periodicDofMap() instead.")]]
188 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
189 { return periodicDofMap(); }
190
192 friend inline LocalView localView(const BoxFVGridGeometry& gg)
193 { return { gg.cache_ }; }
194
195private:
196
197 class BoxGridGeometryCache
198 {
199 friend class BoxFVGridGeometry;
200 public:
202 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
203
204 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
205 : gridGeometry_(&gg)
206 {}
207
208 const BoxFVGridGeometry& gridGeometry() const
209 { return *gridGeometry_; }
210
212 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
213 { return scvs_[eIdx]; }
214
216 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
217 { return scvfs_[eIdx]; }
218
220 bool hasBoundaryScvf(GridIndexType eIdx) const
221 { return hasBoundaryScvf_[eIdx]; }
222
224 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx) const
225 { return scvfBoundaryGeometryKeys_.at(eIdx); }
226
227 private:
228 void clear_()
229 {
230 scvs_.clear();
231 scvfs_.clear();
232 hasBoundaryScvf_.clear();
233 scvfBoundaryGeometryKeys_.clear();
234 }
235
236 std::vector<std::vector<SubControlVolume>> scvs_;
237 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
238 std::vector<bool> hasBoundaryScvf_;
239 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
240
241 const BoxFVGridGeometry* gridGeometry_;
242 };
243
244public:
247 using Cache = BoxGridGeometryCache;
248
249private:
250 using GeometryHelper = typename Cache::GeometryHelper;
251
252 void update_()
253 {
254 cache_.clear_();
255
256 const auto numElements = this->gridView().size(0);
257 cache_.scvs_.resize(numElements);
258 cache_.scvfs_.resize(numElements);
259 cache_.hasBoundaryScvf_.resize(numElements, false);
260
261 boundaryDofIndices_.assign(numDofs(), false);
262
263 numScv_ = 0;
264 numScvf_ = 0;
265 numBoundaryScvf_ = 0;
266 // Build the SCV and SCV faces
267 for (const auto& element : elements(this->gridView()))
268 {
269 // fill the element map with seeds
270 const auto eIdx = this->elementMapper().index(element);
271
272 // count
273 numScv_ += element.subEntities(dim);
274 numScvf_ += element.subEntities(dim-1);
275
276 // get the element geometry
277 auto elementGeometry = element.geometry();
278 const auto refElement = referenceElement(elementGeometry);
279
280 // instantiate the geometry helper
281 GeometryHelper geometryHelper(elementGeometry);
282
283 // construct the sub control volumes
284 cache_.scvs_[eIdx].resize(elementGeometry.corners());
285 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
286 {
287 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
288
289 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(
290 geometryHelper.getScvCorners(scvLocalIdx),
291 scvLocalIdx,
292 eIdx,
293 dofIdxGlobal
294 );
295 }
296
297 // construct the sub control volume faces
298 LocalIndexType scvfLocalIdx = 0;
299 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
300 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
301 {
302 // find the global and local scv indices this scvf is belonging to
303 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
304 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
305
306 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
307 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(
308 corners,
309 geometryHelper.normal(corners, localScvIndices),
310 element,
311 elementGeometry,
312 scvfLocalIdx,
313 std::move(localScvIndices),
314 false
315 );
316 }
317
318 // construct the sub control volume faces on the domain boundary
319 for (const auto& intersection : intersections(this->gridView(), element))
320 {
321 if (intersection.boundary() && !intersection.neighbor())
322 {
323 const auto isGeometry = intersection.geometry();
324 cache_.hasBoundaryScvf_[eIdx] = true;
325
326 // count
327 numScvf_ += isGeometry.corners();
328 numBoundaryScvf_ += isGeometry.corners();
329
330 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
331 {
332 // find the scvs this scvf is belonging to
333 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
334 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
335
336 cache_.scvfs_[eIdx].emplace_back(
337 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
338 intersection.centerUnitOuterNormal(),
339 intersection,
340 isGeometry,
341 isScvfLocalIdx,
342 scvfLocalIdx,
343 std::move(localScvIndices),
344 true
345 );
346
347 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
348 static_cast<LocalIndexType>(intersection.indexInInside()),
349 static_cast<LocalIndexType>(isScvfLocalIdx)
350 }});
351
352 // increment local counter
353 scvfLocalIdx++;
354 }
355
356 // add all vertices on the intersection to the set of
357 // boundary vertices
358 const auto fIdx = intersection.indexInInside();
359 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
360 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
361 {
362 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
363 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
364 boundaryDofIndices_[vIdxGlobal] = true;
365 }
366 }
367
368 // inform the grid geometry if we have periodic boundaries
369 else if (intersection.boundary() && intersection.neighbor())
370 {
371 this->setPeriodic();
372
373 // find the mapped periodic vertex of all vertices on periodic boundaries
374 const auto fIdx = intersection.indexInInside();
375 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
376 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
377 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
378 {
379 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
380 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
381 const auto vPos = elementGeometry.corner(vIdx);
382
383 const auto& outside = intersection.outside();
384 const auto outsideGeometry = outside.geometry();
385 for (const auto& isOutside : intersections(this->gridView(), outside))
386 {
387 // only check periodic vertices of the periodic neighbor
388 if (isOutside.boundary() && isOutside.neighbor())
389 {
390 const auto fIdxOutside = isOutside.indexInInside();
391 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
392 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
393 {
394 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
395 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
396 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
397 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
398 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
399 }
400 }
401 }
402 }
403 }
404 }
405 }
406
407 // error check: periodic boundaries currently don't work for box in parallel
408 if (this->isPeriodic() && this->gridView().comm().size() > 1)
409 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
410 }
411
412 const FeCache feCache_;
413
414 std::size_t numScv_;
415 std::size_t numScvf_;
416 std::size_t numBoundaryScvf_;
417
418 // vertices on the boundary
419 std::vector<bool> boundaryDofIndices_;
420
421 // a map for periodic boundary vertices
422 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
423
424 Cache cache_;
425};
426
434template<class Scalar, class GV, class Traits>
435class BoxFVGridGeometry<Scalar, GV, false, Traits>
436: public BaseGridGeometry<GV, Traits>
437{
440 using GridIndexType = typename IndexTraits<GV>::GridIndex;
441
442 static const int dim = GV::dimension;
443 static const int dimWorld = GV::dimensionworld;
444
445 using Element = typename GV::template Codim<0>::Entity;
446 using CoordScalar = typename GV::ctype;
447
448public:
451 static constexpr DiscretizationMethod discMethod{};
452
456 using LocalView = typename Traits::template LocalView<ThisType, false>;
458 using SubControlVolume = typename Traits::SubControlVolume;
460 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
464 using DofMapper = typename Traits::VertexMapper;
466 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
468 using GridView = GV;
469
471 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
472 : ParentType(std::move(gg))
473 , cache_(*this)
474 {
475 update_();
476 }
477
480 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
481 {}
482
485 const DofMapper& dofMapper() const
486 { return this->vertexMapper(); }
487
489 std::size_t numScv() const
490 { return numScv_; }
491
493 std::size_t numScvf() const
494 { return numScvf_; }
495
498 std::size_t numBoundaryScvf() const
499 { return numBoundaryScvf_; }
500
502 std::size_t numDofs() const
503 { return this->vertexMapper().size(); }
504
505
507 void update(const GridView& gridView)
508 {
509 ParentType::update(gridView);
510 update_();
511 }
512
514 void update(GridView&& gridView)
515 {
516 ParentType::update(std::move(gridView));
517 update_();
518 }
519
521 const FeCache& feCache() const
522 { return feCache_; }
523
525 bool dofOnBoundary(GridIndexType dofIdx) const
526 { return boundaryDofIndices_[dofIdx]; }
527
529 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
530 { return periodicVertexMap_.count(dofIdx); }
531
533 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
534 { return periodicVertexMap_.at(dofIdx); }
535
537 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
538 { return periodicVertexMap_; }
539
541 [[deprecated("Will be removed after release 3.9. Use periodicDofMap() instead.")]]
542 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
543 { return periodicDofMap(); }
544
546 friend inline LocalView localView(const BoxFVGridGeometry& gg)
547 { return { gg.cache_ }; }
548
549private:
550
551 class BoxGridGeometryCache
552 {
553 friend class BoxFVGridGeometry;
554 public:
556 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
557
558 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
559 : gridGeometry_(&gg)
560 {}
561
562 const BoxFVGridGeometry& gridGeometry() const
563 { return *gridGeometry_; }
564
565 private:
566 const BoxFVGridGeometry* gridGeometry_;
567 };
568
569public:
572 using Cache = BoxGridGeometryCache;
573
574private:
575
576 void update_()
577 {
578 boundaryDofIndices_.assign(numDofs(), false);
579
580 // save global data on the grid's scvs and scvfs
581 // TODO do we need those information?
582 numScv_ = 0;
583 numScvf_ = 0;
584 numBoundaryScvf_ = 0;
585 for (const auto& element : elements(this->gridView()))
586 {
587 numScv_ += element.subEntities(dim);
588 numScvf_ += element.subEntities(dim-1);
589
590 const auto elementGeometry = element.geometry();
591 const auto refElement = referenceElement(elementGeometry);
592
593 // store the sub control volume face indices on the domain boundary
594 for (const auto& intersection : intersections(this->gridView(), element))
595 {
596 if (intersection.boundary() && !intersection.neighbor())
597 {
598 const auto isGeometry = intersection.geometry();
599 numScvf_ += isGeometry.corners();
600 numBoundaryScvf_ += isGeometry.corners();
601
602 // add all vertices on the intersection to the set of
603 // boundary vertices
604 const auto fIdx = intersection.indexInInside();
605 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
606 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
607 {
608 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
609 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
610 boundaryDofIndices_[vIdxGlobal] = true;
611 }
612 }
613
614 // inform the grid geometry if we have periodic boundaries
615 else if (intersection.boundary() && intersection.neighbor())
616 {
617 this->setPeriodic();
618
619 // find the mapped periodic vertex of all vertices on periodic boundaries
620 const auto fIdx = intersection.indexInInside();
621 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
622 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
623 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
624 {
625 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
626 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
627 const auto vPos = elementGeometry.corner(vIdx);
628
629 const auto& outside = intersection.outside();
630 const auto outsideGeometry = outside.geometry();
631 for (const auto& isOutside : intersections(this->gridView(), outside))
632 {
633 // only check periodic vertices of the periodic neighbor
634 if (isOutside.boundary() && isOutside.neighbor())
635 {
636 const auto fIdxOutside = isOutside.indexInInside();
637 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
638 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
639 {
640 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
641 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
642 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
643 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
644 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
645 }
646 }
647 }
648 }
649 }
650 }
651 }
652
653 // error check: periodic boundaries currently don't work for box in parallel
654 if (this->isPeriodic() && this->gridView().comm().size() > 1)
655 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
656 }
657
658 const FeCache feCache_;
659
660 // Information on the global number of geometries
661 // TODO do we need those information?
662 std::size_t numScv_;
663 std::size_t numScvf_;
664 std::size_t numBoundaryScvf_;
665
666 // vertices on the boundary
667 std::vector<bool> boundaryDofIndices_;
668
669 // a map for periodic boundary vertices
670 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
671
672 Cache cache_;
673};
674
675} // end namespace Dumux
676
677#endif
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:41
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:437
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:498
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:546
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:529
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:456
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:460
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:479
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:458
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:533
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:521
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:572
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:493
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:464
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:454
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:525
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:485
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:514
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:462
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:502
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:542
BoxFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg)
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/box/fvgridgeometry.hh:471
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:489
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:507
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:466
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:537
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:83
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:192
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:247
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:188
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:153
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:179
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:148
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:135
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:100
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:167
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:144
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:112
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:106
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:131
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:183
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:125
BoxFVGridGeometry(std::shared_ptr< BasicGridGeometry > gg)
Constructor with basic grid geometry used to share state with another grid geometry on the same grid ...
Definition: discretization/box/fvgridgeometry.hh:117
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:108
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:102
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:171
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:139
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:175
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:160
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:104
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:110
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:72
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:257
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:62
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:60
Defines the default element and vertex mapper types.
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.
Helper classes to compute the integration elements.
Dune::Std::detected_or_t< Dumux::BasicGridGeometry< GV, typename T::ElementMapper, typename T::VertexMapper >, Detail::SpecifiesBaseGridGeometry, T > BasicGridGeometry_t
Type of the basic grid geometry implementation used as backend.
Definition: basegridgeometry.hh:42
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Dune::Std::detected_or_t< Dumux::BoxGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxGeometryHelper_t
Definition: discretization/box/fvgridgeometry.hh:42
typename T::GeometryHelper SpecifiesGeometryHelper
Definition: basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition: discretization/box/fvgridgeometry.hh:54
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26