version 3.11-dev
Loading...
Searching...
No Matches
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#include <span>
22
23#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
24
35
37
38namespace Dumux {
39
40namespace Detail {
41template<class GV, class T>
42using BoxGeometryHelper_t = Dune::Std::detected_or_t<
45 T
46>;
47} // end namespace Detail
48
53template<class GridView, class ScvRule = Dumux::QuadratureRules::MidpointQuadrature, class ScvfRule = Dumux::QuadratureRules::MidpointQuadrature>
55
62template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>, class QuadratureTraits = BoxQuadratureTraits<GridView>>
64: public MapperTraits, public QuadratureTraits
65{
68
69 template<class GridGeometry, bool enableCache>
71};
72
79template<class Scalar,
80 class GridView,
81 bool enableGridGeometryCache = false,
84
91template<class Scalar, class GV, class Traits>
92class BoxFVGridGeometry<Scalar, GV, true, Traits>
93: public BaseGridGeometry<GV, Traits>
94{
96 using ParentType = BaseGridGeometry<GV, Traits>;
97 using GridIndexType = typename IndexTraits<GV>::GridIndex;
98 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
99
100 using Element = typename GV::template Codim<0>::Entity;
101 using CoordScalar = typename GV::ctype;
102 static const int dim = GV::dimension;
103 static const int dimWorld = GV::dimensionworld;
104
105public:
109
113 using LocalView = typename Traits::template LocalView<ThisType, true>;
115 using SubControlVolume = typename Traits::SubControlVolume;
117 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
123 using DofMapper = typename Traits::VertexMapper;
125 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
127 using GridView = GV;
131 using ScvQuadratureRule = typename Traits::ScvQuadratureRule;
133 using ScvfQuadratureRule = typename Traits::ScvfQuadratureRule;
134
136 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
137 : ParentType(std::move(gg))
138 , cache_(*this)
139 , periodicGridTraits_(this->gridView().grid())
140 {
141 update_();
142 }
143
148
151 const DofMapper& dofMapper() const
152 { return this->vertexMapper(); }
153
155 std::size_t numScv() const
156 { return numScv_; }
157
159 std::size_t numScvf() const
160 { return numScvf_; }
161
164 std::size_t numBoundaryScvf() const
165 { return numBoundaryScvf_; }
166
168 std::size_t numDofs() const
169 { return this->vertexMapper().size(); }
170
171
174 {
176 update_();
177 }
178
181 {
182 ParentType::update(std::move(gridView));
183 update_();
184 }
185
187 const FeCache& feCache() const
188 { return feCache_; }
189
191 bool dofOnBoundary(GridIndexType dofIdx) const
192 { return boundaryDofIndices_[dofIdx]; }
193
195 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
196 { return periodicDofMap_.count(dofIdx); }
197
199 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
200 { return periodicDofMap_.at(dofIdx); }
201
203 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
204 { return periodicDofMap_; }
205
207 friend inline LocalView localView(const BoxFVGridGeometry& gg)
208 { return { gg.cache_ }; }
209
210private:
211
212 class BoxGridGeometryCache
213 {
214 friend class BoxFVGridGeometry;
215 public:
217 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
218
219 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
220 : gridGeometry_(&gg)
221 {}
222
223 const BoxFVGridGeometry& gridGeometry() const
224 { return *gridGeometry_; }
225
227 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
228 { return scvs_[eIdx]; }
229
231 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
232 { return scvfs_[eIdx]; }
233
235 bool hasBoundaryScvf(GridIndexType eIdx) const
236 { return hasBoundaryScvf_[eIdx]; }
237
239 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx) const
240 { return scvfBoundaryGeometryKeys_.at(eIdx); }
241
243 auto boundaryFaces(GridIndexType eIdx) const -> std::span<const BoundaryFace>
244 {
245 if (auto it = boundaryFaces_.find(eIdx); it != boundaryFaces_.end())
246 return {it->second};
247 return {};
248 }
249
251 const auto& boundaryFaceScvfRanges(GridIndexType eIdx) const
252 { return boundaryFaceScvfRanges_.at(eIdx); }
253
254 private:
255 void clear_()
256 {
257 scvs_.clear();
258 scvfs_.clear();
259 hasBoundaryScvf_.clear();
260 scvfBoundaryGeometryKeys_.clear();
261 boundaryFaces_.clear();
262 boundaryFaceScvfRanges_.clear();
263 }
264
265 std::vector<std::vector<SubControlVolume>> scvs_;
266 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
267 std::vector<bool> hasBoundaryScvf_;
268 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
269 std::unordered_map<GridIndexType, Dune::ReservedVector<typename BoxFVGridGeometry::BoundaryFace, 2*dim>> boundaryFaces_;
270 std::unordered_map<GridIndexType, Dune::ReservedVector<std::array<LocalIndexType, 2>, 2*dim>> boundaryFaceScvfRanges_;
271
272 const BoxFVGridGeometry* gridGeometry_;
273 };
274
275public:
278 using Cache = BoxGridGeometryCache;
279
280private:
281 using GeometryHelper = typename Cache::GeometryHelper;
282
283 void update_()
284 {
285 cache_.clear_();
286
287 const auto numElements = this->gridView().size(0);
288 cache_.scvs_.resize(numElements);
289 cache_.scvfs_.resize(numElements);
290 cache_.hasBoundaryScvf_.resize(numElements, false);
291
292 boundaryDofIndices_.assign(numDofs(), false);
293
294 numScv_ = 0;
295 numScvf_ = 0;
296 numBoundaryScvf_ = 0;
297 // Build the SCV and SCV faces
298 for (const auto& element : elements(this->gridView()))
299 {
300 // fill the element map with seeds
301 const auto eIdx = this->elementMapper().index(element);
302
303 // count
304 numScv_ += element.subEntities(dim);
305 numScvf_ += element.subEntities(dim-1);
306
307 // get the element geometry
308 auto elementGeometry = element.geometry();
309 const auto refElement = referenceElement(elementGeometry);
310
311 // instantiate the geometry helper
312 GeometryHelper geometryHelper(elementGeometry);
313
314 // construct the sub control volumes
315 cache_.scvs_[eIdx].resize(elementGeometry.corners());
316 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
317 {
318 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
319
320 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(
321 geometryHelper.getScvCorners(scvLocalIdx),
322 scvLocalIdx,
323 eIdx,
324 dofIdxGlobal
325 );
326 }
327
328 // construct the sub control volume faces
329 LocalIndexType scvfLocalIdx = 0;
330 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
331 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
332 {
333 // find the global and local scv indices this scvf is belonging to
334 std::array<LocalIndexType, 2> localScvIndices{{
335 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
336 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
337 }};
338
339 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
340 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(
341 corners,
342 geometryHelper.normal(corners, localScvIndices),
343 element,
344 scvfLocalIdx,
345 std::move(localScvIndices)
346 );
347 }
348
349 // construct the sub control volume faces on the domain boundary
350 LocalIndexType numBoundaryFaces = 0;
351 for (const auto& intersection : intersections(this->gridView(), element))
352 {
353 if (intersection.boundary() && !intersection.neighbor())
354 {
355 const auto isGeometry = intersection.geometry();
356 cache_.hasBoundaryScvf_[eIdx] = true;
357
358 // add one boundary face per boundary intersection
359 cache_.boundaryFaces_[eIdx].push_back(BoundaryFace{
360 isGeometry.center(),
361 isGeometry.volume(),
362 intersection.centerUnitOuterNormal(),
363 numBoundaryFaces++,
364 static_cast<LocalIndexType>(intersection.indexInInside()),
365 typename BoundaryFace::Traits::BoundaryFlag{intersection}
366 });
367
368 // count
369 numScvf_ += isGeometry.corners();
370 numBoundaryScvf_ += isGeometry.corners();
371
372 // record the scvf subrange for this boundary face: {offset, count}
373 cache_.boundaryFaceScvfRanges_[eIdx].push_back(std::array<LocalIndexType, 2>{{
374 scvfLocalIdx,
375 static_cast<LocalIndexType>(isGeometry.corners())
376 }});
377
378 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
379 {
380 // find the scvs this scvf is belonging to
381 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
382 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
383
384 cache_.scvfs_[eIdx].emplace_back(
385 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
386 intersection.centerUnitOuterNormal(),
387 intersection,
388 isScvfLocalIdx,
389 scvfLocalIdx,
390 std::move(localScvIndices)
391 );
392
393 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
394 static_cast<LocalIndexType>(intersection.indexInInside()),
395 static_cast<LocalIndexType>(isScvfLocalIdx)
396 }});
397
398 // increment local counter
399 scvfLocalIdx++;
400 }
401
402 // add all vertices on the intersection to the set of
403 // boundary vertices
404 const auto fIdx = intersection.indexInInside();
405 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
406 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
407 {
408 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
409 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
410 boundaryDofIndices_[vIdxGlobal] = true;
411 }
412 }
413
414 // inform the grid geometry if we have periodic boundaries
415 else if (periodicGridTraits_.isPeriodic(intersection))
416 {
417 this->setPeriodic();
418
419 // find the mapped periodic vertex of all vertices on periodic boundaries
420 const auto fIdx = intersection.indexInInside();
421 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
422 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
423 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
424 {
425 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
426 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
427 const auto vPos = elementGeometry.corner(vIdx);
428
429 const auto& outside = intersection.outside();
430 const auto outsideGeometry = outside.geometry();
431 for (const auto& isOutside : intersections(this->gridView(), outside))
432 {
433 // only check periodic vertices of the periodic neighbor
434 if (periodicGridTraits_.isPeriodic(isOutside))
435 {
436 const auto fIdxOutside = isOutside.indexInInside();
437 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
438 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
439 {
440 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
441 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
442 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
443 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
444 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
445 }
446 }
447 }
448 }
449 }
450 }
451 }
452
453 // error check: periodic boundaries currently don't work for box in parallel
454 if (this->isPeriodic() && this->gridView().comm().size() > 1)
455 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
456 }
457
458 const FeCache feCache_;
459
460 std::size_t numScv_;
461 std::size_t numScvf_;
462 std::size_t numBoundaryScvf_;
463
464 // vertices on the boundary
465 std::vector<bool> boundaryDofIndices_;
466
467 // a map for periodic boundary vertices
468 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
469
470 Cache cache_;
471
473};
474
482template<class Scalar, class GV, class Traits>
483class BoxFVGridGeometry<Scalar, GV, false, Traits>
484: public BaseGridGeometry<GV, Traits>
485{
487 using ParentType = BaseGridGeometry<GV, Traits>;
488 using GridIndexType = typename IndexTraits<GV>::GridIndex;
489
490 static const int dim = GV::dimension;
491 static const int dimWorld = GV::dimensionworld;
492
493 using Element = typename GV::template Codim<0>::Entity;
494 using CoordScalar = typename GV::ctype;
495
496public:
500
504 using LocalView = typename Traits::template LocalView<ThisType, false>;
506 using SubControlVolume = typename Traits::SubControlVolume;
508 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
514 using DofMapper = typename Traits::VertexMapper;
516 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
518 using GridView = GV;
522 using ScvQuadratureRule = typename Traits::ScvQuadratureRule;
524 using ScvfQuadratureRule = typename Traits::ScvfQuadratureRule;
525
527 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
528 : ParentType(std::move(gg))
529 , cache_(*this)
530 , periodicGridTraits_(this->gridView().grid())
531 {
532 update_();
533 }
534
539
542 const DofMapper& dofMapper() const
543 { return this->vertexMapper(); }
544
546 std::size_t numScv() const
547 { return numScv_; }
548
550 std::size_t numScvf() const
551 { return numScvf_; }
552
555 std::size_t numBoundaryScvf() const
556 { return numBoundaryScvf_; }
557
559 std::size_t numDofs() const
560 { return this->vertexMapper().size(); }
561
562
565 {
567 update_();
568 }
569
572 {
573 ParentType::update(std::move(gridView));
574 update_();
575 }
576
578 const FeCache& feCache() const
579 { return feCache_; }
580
582 bool dofOnBoundary(GridIndexType dofIdx) const
583 { return boundaryDofIndices_[dofIdx]; }
584
586 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
587 { return periodicDofMap_.count(dofIdx); }
588
590 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
591 { return periodicDofMap_.at(dofIdx); }
592
594 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
595 { return periodicDofMap_; }
596
598 friend inline LocalView localView(const BoxFVGridGeometry& gg)
599 { return { gg.cache_ }; }
600
601private:
602
603 class BoxGridGeometryCache
604 {
605 friend class BoxFVGridGeometry;
606 public:
608 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
609
610 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
611 : gridGeometry_(&gg)
612 {}
613
614 const BoxFVGridGeometry& gridGeometry() const
615 { return *gridGeometry_; }
616
617 private:
618 const BoxFVGridGeometry* gridGeometry_;
619 };
620
621public:
624 using Cache = BoxGridGeometryCache;
625
626private:
627
628 void update_()
629 {
630 boundaryDofIndices_.assign(numDofs(), false);
631
632 // save global data on the grid's scvs and scvfs
633 // TODO do we need those information?
634 numScv_ = 0;
635 numScvf_ = 0;
636 numBoundaryScvf_ = 0;
637 for (const auto& element : elements(this->gridView()))
638 {
639 numScv_ += element.subEntities(dim);
640 numScvf_ += element.subEntities(dim-1);
641
642 const auto elementGeometry = element.geometry();
643 const auto refElement = referenceElement(elementGeometry);
644
645 // store the sub control volume face indices on the domain boundary
646 for (const auto& intersection : intersections(this->gridView(), element))
647 {
648 if (intersection.boundary() && !intersection.neighbor())
649 {
650 const auto isGeometry = intersection.geometry();
651 numScvf_ += isGeometry.corners();
652 numBoundaryScvf_ += isGeometry.corners();
653
654 // add all vertices on the intersection to the set of
655 // boundary vertices
656 const auto fIdx = intersection.indexInInside();
657 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
658 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
659 {
660 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
661 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
662 boundaryDofIndices_[vIdxGlobal] = true;
663 }
664 }
665
666 // inform the grid geometry if we have periodic boundaries
667 else if (periodicGridTraits_.isPeriodic(intersection))
668 {
669 this->setPeriodic();
670
671 // find the mapped periodic vertex of all vertices on periodic boundaries
672 const auto fIdx = intersection.indexInInside();
673 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
674 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
675 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
676 {
677 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
678 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
679 const auto vPos = elementGeometry.corner(vIdx);
680
681 const auto& outside = intersection.outside();
682 const auto outsideGeometry = outside.geometry();
683 for (const auto& isOutside : intersections(this->gridView(), outside))
684 {
685 // only check periodic vertices of the periodic neighbor
686 if (periodicGridTraits_.isPeriodic(isOutside))
687 {
688 const auto fIdxOutside = isOutside.indexInInside();
689 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
690 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
691 {
692 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
693 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
694 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
695 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
696 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
697 }
698 }
699 }
700 }
701 }
702 }
703 }
704
705 // error check: periodic boundaries currently don't work for box in parallel
706 if (this->isPeriodic() && this->gridView().comm().size() > 1)
707 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
708 }
709
710 const FeCache feCache_;
711
712 // Information on the global number of geometries
713 // TODO do we need those information?
714 std::size_t numScv_;
715 std::size_t numScvf_;
716 std::size_t numBoundaryScvf_;
717
718 // vertices on the boundary
719 std::vector<bool> boundaryDofIndices_;
720
721 // a map for periodic boundary vertices
722 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
723
724 Cache cache_;
725
727};
728
729} // end namespace Dumux
730
731#endif
Base class for grid geometries.
Implementation of a boundary face related to primary grid elements (dune intersections).
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices for constant grids.
Definition basegridgeometry.hh:112
void setPeriodic(bool value=true)
Set the periodicity of the grid geometry.
Definition basegridgeometry.hh:169
const GlobalCoordinate & bBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition basegridgeometry.hh:156
Element element(GridIndexType eIdx) const
Get an element from a global element index.
Definition basegridgeometry.hh:142
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices for constant grids.
Definition basegridgeometry.hh:106
const GridView & gridView() const
Return the gridView this grid geometry object lives on.
Definition basegridgeometry.hh:100
void update(const GridView &gridView)
Update all fvElementGeometries (call this after grid adaption).
Definition basegridgeometry.hh:88
const GlobalCoordinate & bBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition basegridgeometry.hh:149
bool isPeriodic() const
Returns if the grid geometry is periodic (at all).
Definition basegridgeometry.hh:162
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition discretization/box/fvelementgeometry.hh:46
std::size_t numBoundaryScvf() const
Definition discretization/box/fvgridgeometry.hh:555
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition discretization/box/fvgridgeometry.hh:598
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition discretization/box/fvgridgeometry.hh:586
GV GridView
export the grid view type
Definition discretization/box/fvgridgeometry.hh:518
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition discretization/box/fvgridgeometry.hh:504
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:508
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition discretization/box/fvgridgeometry.hh:536
static constexpr DiscretizationMethod discMethod
Definition discretization/box/fvgridgeometry.hh:499
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:506
Experimental::BoundaryFace< GV > BoundaryFace
export the boundary face type
Definition discretization/box/fvgridgeometry.hh:510
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvgridgeometry.hh:522
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:590
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition discretization/box/fvgridgeometry.hh:578
BoxGridGeometryCache Cache
Definition discretization/box/fvgridgeometry.hh:624
std::size_t numScvf() const
The total number of sun control volume faces.
Definition discretization/box/fvgridgeometry.hh:550
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition discretization/box/fvgridgeometry.hh:514
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition discretization/box/fvgridgeometry.hh:502
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition discretization/box/fvgridgeometry.hh:582
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvgridgeometry.hh:524
const DofMapper & dofMapper() const
Definition discretization/box/fvgridgeometry.hh:542
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition discretization/box/fvgridgeometry.hh:498
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition discretization/box/fvgridgeometry.hh:520
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:571
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/box/fvgridgeometry.hh:512
std::size_t numDofs() const
The total number of degrees of freedom.
Definition discretization/box/fvgridgeometry.hh:559
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:527
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvgridgeometry.hh:546
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:564
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition discretization/box/fvgridgeometry.hh:516
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition discretization/box/fvgridgeometry.hh:594
GV GridView
export the grid view type
Definition discretization/box/fvgridgeometry.hh:127
Experimental::BoundaryFace< GV > BoundaryFace
export the boundary face type
Definition discretization/box/fvgridgeometry.hh:119
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition discretization/box/fvgridgeometry.hh:207
BoxGridGeometryCache Cache
Definition discretization/box/fvgridgeometry.hh:278
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:173
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:199
std::size_t numDofs() const
The total number of degrees of freedom.
Definition discretization/box/fvgridgeometry.hh:168
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvgridgeometry.hh:155
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition discretization/box/fvgridgeometry.hh:111
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvgridgeometry.hh:131
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition discretization/box/fvgridgeometry.hh:187
std::size_t numBoundaryScvf() const
Definition discretization/box/fvgridgeometry.hh:164
DiscretizationMethods::Box DiscretizationMethod
export the discretization method this geometry belongs to
Definition discretization/box/fvgridgeometry.hh:107
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvgridgeometry.hh:133
static constexpr DiscretizationMethod discMethod
Definition discretization/box/fvgridgeometry.hh:108
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition discretization/box/fvgridgeometry.hh:125
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:117
const DofMapper & dofMapper() const
Definition discretization/box/fvgridgeometry.hh:151
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition discretization/box/fvgridgeometry.hh:203
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition discretization/box/fvgridgeometry.hh:145
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:136
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition discretization/box/fvgridgeometry.hh:121
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition discretization/box/fvgridgeometry.hh:129
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition discretization/box/fvgridgeometry.hh:113
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition discretization/box/fvgridgeometry.hh:191
std::size_t numScvf() const
The total number of sun control volume faces.
Definition discretization/box/fvgridgeometry.hh:159
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition discretization/box/fvgridgeometry.hh:195
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition discretization/box/fvgridgeometry.hh:180
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition discretization/box/fvgridgeometry.hh:115
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition discretization/box/fvgridgeometry.hh:123
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition discretization/box/fvgridgeometry.hh:83
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
Class for a boundary face related to primary grid elements (dune intersections).
Definition boundaryface.hh:69
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.
CVFE::DefaultQuadratureTraits< GridView, ScvRule, ScvfRule > BoxQuadratureTraits
Quadrature rule traits for Box discretization.
Definition discretization/box/fvgridgeometry.hh:54
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:38
BaseGridGeometry(std::shared_ptr< BaseImplementation > impl)
Constructor from a BaseImplementation.
Definition basegridgeometry.hh:72
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Definition cvfelocalresidual.hh:25
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:102
Definition adapt.hh:17
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
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:65
BoxFVElementGeometry< GridGeometry, enableCache > LocalView
Definition discretization/box/fvgridgeometry.hh:70
BoxSubControlVolume< GridView > SubControlVolume
Definition discretization/box/fvgridgeometry.hh:66
BoxSubControlVolumeFace< GridView > SubControlVolumeFace
Definition discretization/box/fvgridgeometry.hh:67
Quadrature rule traits for discretization schemes.
Definition quadraturerules.hh:85
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:27
unsigned int LocalIndex
Definition indextraits.hh:28
Definition periodicgridtraits.hh:24
Definition periodicgridtraits.hh:23