version 3.11-dev
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-FileCopyrightText: 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
35
36namespace Dumux {
37
38namespace Detail {
39template<class GV, class T>
40using BoxGeometryHelper_t = Dune::Std::detected_or_t<
43 T
44>;
45} // end namespace Detail
46
53template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
55: public MapperTraits
56{
59
60 template<class GridGeometry, bool enableCache>
62};
63
70template<class Scalar,
71 class GridView,
72 bool enableGridGeometryCache = false,
75
82template<class Scalar, class GV, class Traits>
83class BoxFVGridGeometry<Scalar, GV, true, Traits>
84: public BaseGridGeometry<GV, Traits>
85{
88 using GridIndexType = typename IndexTraits<GV>::GridIndex;
89 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
90
91 using Element = typename GV::template Codim<0>::Entity;
92 using CoordScalar = typename GV::ctype;
93 static const int dim = GV::dimension;
94 static const int dimWorld = GV::dimensionworld;
95
96public:
99 static constexpr DiscretizationMethod discMethod{};
100
104 using LocalView = typename Traits::template LocalView<ThisType, true>;
106 using SubControlVolume = typename Traits::SubControlVolume;
108 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
112 using DofMapper = typename Traits::VertexMapper;
114 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
116 using GridView = GV;
119
121 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
122 : ParentType(std::move(gg))
123 , cache_(*this)
124 , periodicGridTraits_(this->gridView().grid())
125 {
126 update_();
127 }
128
131 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
132 {}
133
136 const DofMapper& dofMapper() const
137 { return this->vertexMapper(); }
138
140 std::size_t numScv() const
141 { return numScv_; }
142
144 std::size_t numScvf() const
145 { return numScvf_; }
146
149 std::size_t numBoundaryScvf() const
150 { return numBoundaryScvf_; }
151
153 std::size_t numDofs() const
154 { return this->vertexMapper().size(); }
155
156
158 void update(const GridView& gridView)
159 {
160 ParentType::update(gridView);
161 update_();
162 }
163
165 void update(GridView&& gridView)
166 {
167 ParentType::update(std::move(gridView));
168 update_();
169 }
170
172 const FeCache& feCache() const
173 { return feCache_; }
174
176 bool dofOnBoundary(GridIndexType dofIdx) const
177 { return boundaryDofIndices_[dofIdx]; }
178
180 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
181 { return periodicDofMap_.count(dofIdx); }
182
184 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
185 { return periodicDofMap_.at(dofIdx); }
186
188 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() 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::array<LocalIndexType, 2> localScvIndices{{
304 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
305 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
306 }};
307
308 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
309 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(
310 corners,
311 geometryHelper.normal(corners, localScvIndices),
312 element,
313 scvfLocalIdx,
314 std::move(localScvIndices)
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::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
335
336 cache_.scvfs_[eIdx].emplace_back(
337 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
338 intersection.centerUnitOuterNormal(),
339 intersection,
340 isScvfLocalIdx,
341 scvfLocalIdx,
342 std::move(localScvIndices)
343 );
344
345 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
346 static_cast<LocalIndexType>(intersection.indexInInside()),
347 static_cast<LocalIndexType>(isScvfLocalIdx)
348 }});
349
350 // increment local counter
351 scvfLocalIdx++;
352 }
353
354 // add all vertices on the intersection to the set of
355 // boundary vertices
356 const auto fIdx = intersection.indexInInside();
357 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
358 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
359 {
360 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
361 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
362 boundaryDofIndices_[vIdxGlobal] = true;
363 }
364 }
365
366 // inform the grid geometry if we have periodic boundaries
367 else if (periodicGridTraits_.isPeriodic(intersection))
368 {
369 this->setPeriodic();
370
371 // find the mapped periodic vertex of all vertices on periodic boundaries
372 const auto fIdx = intersection.indexInInside();
373 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
374 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
375 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
376 {
377 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
378 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
379 const auto vPos = elementGeometry.corner(vIdx);
380
381 const auto& outside = intersection.outside();
382 const auto outsideGeometry = outside.geometry();
383 for (const auto& isOutside : intersections(this->gridView(), outside))
384 {
385 // only check periodic vertices of the periodic neighbor
386 if (periodicGridTraits_.isPeriodic(isOutside))
387 {
388 const auto fIdxOutside = isOutside.indexInInside();
389 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
390 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
391 {
392 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
393 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
394 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
395 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
396 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
397 }
398 }
399 }
400 }
401 }
402 }
403 }
404
405 // error check: periodic boundaries currently don't work for box in parallel
406 if (this->isPeriodic() && this->gridView().comm().size() > 1)
407 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
408 }
409
410 const FeCache feCache_;
411
412 std::size_t numScv_;
413 std::size_t numScvf_;
414 std::size_t numBoundaryScvf_;
415
416 // vertices on the boundary
417 std::vector<bool> boundaryDofIndices_;
418
419 // a map for periodic boundary vertices
420 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
421
422 Cache cache_;
423
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;
471
473 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
474 : ParentType(std::move(gg))
475 , cache_(*this)
476 , periodicGridTraits_(this->gridView().grid())
477 {
478 update_();
479 }
480
483 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
484 {}
485
488 const DofMapper& dofMapper() const
489 { return this->vertexMapper(); }
490
492 std::size_t numScv() const
493 { return numScv_; }
494
496 std::size_t numScvf() const
497 { return numScvf_; }
498
501 std::size_t numBoundaryScvf() const
502 { return numBoundaryScvf_; }
503
505 std::size_t numDofs() const
506 { return this->vertexMapper().size(); }
507
508
510 void update(const GridView& gridView)
511 {
512 ParentType::update(gridView);
513 update_();
514 }
515
517 void update(GridView&& gridView)
518 {
519 ParentType::update(std::move(gridView));
520 update_();
521 }
522
524 const FeCache& feCache() const
525 { return feCache_; }
526
528 bool dofOnBoundary(GridIndexType dofIdx) const
529 { return boundaryDofIndices_[dofIdx]; }
530
532 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
533 { return periodicDofMap_.count(dofIdx); }
534
536 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
537 { return periodicDofMap_.at(dofIdx); }
538
540 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
541 { return periodicDofMap_; }
542
544 friend inline LocalView localView(const BoxFVGridGeometry& gg)
545 { return { gg.cache_ }; }
546
547private:
548
549 class BoxGridGeometryCache
550 {
551 friend class BoxFVGridGeometry;
552 public:
554 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
555
556 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
557 : gridGeometry_(&gg)
558 {}
559
560 const BoxFVGridGeometry& gridGeometry() const
561 { return *gridGeometry_; }
562
563 private:
564 const BoxFVGridGeometry* gridGeometry_;
565 };
566
567public:
570 using Cache = BoxGridGeometryCache;
571
572private:
573
574 void update_()
575 {
576 boundaryDofIndices_.assign(numDofs(), false);
577
578 // save global data on the grid's scvs and scvfs
579 // TODO do we need those information?
580 numScv_ = 0;
581 numScvf_ = 0;
582 numBoundaryScvf_ = 0;
583 for (const auto& element : elements(this->gridView()))
584 {
585 numScv_ += element.subEntities(dim);
586 numScvf_ += element.subEntities(dim-1);
587
588 const auto elementGeometry = element.geometry();
589 const auto refElement = referenceElement(elementGeometry);
590
591 // store the sub control volume face indices on the domain boundary
592 for (const auto& intersection : intersections(this->gridView(), element))
593 {
594 if (intersection.boundary() && !intersection.neighbor())
595 {
596 const auto isGeometry = intersection.geometry();
597 numScvf_ += isGeometry.corners();
598 numBoundaryScvf_ += isGeometry.corners();
599
600 // add all vertices on the intersection to the set of
601 // boundary vertices
602 const auto fIdx = intersection.indexInInside();
603 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
604 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
605 {
606 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
607 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
608 boundaryDofIndices_[vIdxGlobal] = true;
609 }
610 }
611
612 // inform the grid geometry if we have periodic boundaries
613 else if (periodicGridTraits_.isPeriodic(intersection))
614 {
615 this->setPeriodic();
616
617 // find the mapped periodic vertex of all vertices on periodic boundaries
618 const auto fIdx = intersection.indexInInside();
619 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
620 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
621 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
622 {
623 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
624 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
625 const auto vPos = elementGeometry.corner(vIdx);
626
627 const auto& outside = intersection.outside();
628 const auto outsideGeometry = outside.geometry();
629 for (const auto& isOutside : intersections(this->gridView(), outside))
630 {
631 // only check periodic vertices of the periodic neighbor
632 if (periodicGridTraits_.isPeriodic(isOutside))
633 {
634 const auto fIdxOutside = isOutside.indexInInside();
635 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
636 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
637 {
638 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
639 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
640 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
641 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
642 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
643 }
644 }
645 }
646 }
647 }
648 }
649 }
650
651 // error check: periodic boundaries currently don't work for box in parallel
652 if (this->isPeriodic() && this->gridView().comm().size() > 1)
653 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
654 }
655
656 const FeCache feCache_;
657
658 // Information on the global number of geometries
659 // TODO do we need those information?
660 std::size_t numScv_;
661 std::size_t numScvf_;
662 std::size_t numBoundaryScvf_;
663
664 // vertices on the boundary
665 std::vector<bool> boundaryDofIndices_;
666
667 // a map for periodic boundary vertices
668 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
669
670 Cache cache_;
671
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:42
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:501
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:544
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:532
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:482
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:536
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:524
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:570
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:496
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:528
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:488
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition: discretization/box/fvgridgeometry.hh:470
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:517
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:505
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:473
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:492
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:510
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:540
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:85
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
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:158
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:184
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:153
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:140
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:102
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:172
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:149
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:114
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:108
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:136
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:188
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:130
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:121
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:110
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition: discretization/box/fvgridgeometry.hh:118
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:104
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:176
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:144
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:180
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:165
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:106
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:112
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:74
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:60
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:58
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:44
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
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
Grid properties related to periodicity.
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition: discretization/box/fvgridgeometry.hh:56
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
Definition: periodicgridtraits.hh:24