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
51template<class GridView, class ScvRule = Dumux::QuadratureRules::MidpointQuadrature, class ScvfRule = Dumux::QuadratureRules::MidpointQuadrature>
53
60template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>, class QuadratureTraits = BoxQuadratureTraits<GridView>>
62: public MapperTraits, public QuadratureTraits
63{
66
67 template<class GridGeometry, bool enableCache>
69};
70
77template<class Scalar,
78 class GridView,
79 bool enableGridGeometryCache = false,
82
89template<class Scalar, class GV, class Traits>
90class BoxFVGridGeometry<Scalar, GV, true, Traits>
91: public BaseGridGeometry<GV, Traits>
92{
95 using GridIndexType = typename IndexTraits<GV>::GridIndex;
96 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
97
98 using Element = typename GV::template Codim<0>::Entity;
99 using CoordScalar = typename GV::ctype;
100 static const int dim = GV::dimension;
101 static const int dimWorld = GV::dimensionworld;
102
103public:
106 static constexpr DiscretizationMethod discMethod{};
107
111 using LocalView = typename Traits::template LocalView<ThisType, true>;
113 using SubControlVolume = typename Traits::SubControlVolume;
115 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
119 using DofMapper = typename Traits::VertexMapper;
121 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
123 using GridView = GV;
127 using ScvQuadratureRule = typename Traits::ScvQuadratureRule;
129 using ScvfQuadratureRule = typename Traits::ScvfQuadratureRule;
130
132 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
133 : ParentType(std::move(gg))
134 , cache_(*this)
135 , periodicGridTraits_(this->gridView().grid())
136 {
137 update_();
138 }
139
142 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
143 {}
144
147 const DofMapper& dofMapper() const
148 { return this->vertexMapper(); }
149
151 std::size_t numScv() const
152 { return numScv_; }
153
155 std::size_t numScvf() const
156 { return numScvf_; }
157
160 std::size_t numBoundaryScvf() const
161 { return numBoundaryScvf_; }
162
164 std::size_t numDofs() const
165 { return this->vertexMapper().size(); }
166
167
169 void update(const GridView& gridView)
170 {
171 ParentType::update(gridView);
172 update_();
173 }
174
176 void update(GridView&& gridView)
177 {
178 ParentType::update(std::move(gridView));
179 update_();
180 }
181
183 const FeCache& feCache() const
184 { return feCache_; }
185
187 bool dofOnBoundary(GridIndexType dofIdx) const
188 { return boundaryDofIndices_[dofIdx]; }
189
191 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
192 { return periodicDofMap_.count(dofIdx); }
193
195 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
196 { return periodicDofMap_.at(dofIdx); }
197
199 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
200 { return periodicDofMap_; }
201
203 friend inline LocalView localView(const BoxFVGridGeometry& gg)
204 { return { gg.cache_ }; }
205
206private:
207
208 class BoxGridGeometryCache
209 {
210 friend class BoxFVGridGeometry;
211 public:
213 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
214
215 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
216 : gridGeometry_(&gg)
217 {}
218
219 const BoxFVGridGeometry& gridGeometry() const
220 { return *gridGeometry_; }
221
223 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
224 { return scvs_[eIdx]; }
225
227 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
228 { return scvfs_[eIdx]; }
229
231 bool hasBoundaryScvf(GridIndexType eIdx) const
232 { return hasBoundaryScvf_[eIdx]; }
233
235 const std::vector<std::array<LocalIndexType, 2>>& scvfBoundaryGeometryKeys(GridIndexType eIdx) const
236 { return scvfBoundaryGeometryKeys_.at(eIdx); }
237
238 private:
239 void clear_()
240 {
241 scvs_.clear();
242 scvfs_.clear();
243 hasBoundaryScvf_.clear();
244 scvfBoundaryGeometryKeys_.clear();
245 }
246
247 std::vector<std::vector<SubControlVolume>> scvs_;
248 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
249 std::vector<bool> hasBoundaryScvf_;
250 std::unordered_map<GridIndexType, std::vector<std::array<LocalIndexType, 2>>> scvfBoundaryGeometryKeys_;
251
252 const BoxFVGridGeometry* gridGeometry_;
253 };
254
255public:
258 using Cache = BoxGridGeometryCache;
259
260private:
261 using GeometryHelper = typename Cache::GeometryHelper;
262
263 void update_()
264 {
265 cache_.clear_();
266
267 const auto numElements = this->gridView().size(0);
268 cache_.scvs_.resize(numElements);
269 cache_.scvfs_.resize(numElements);
270 cache_.hasBoundaryScvf_.resize(numElements, false);
271
272 boundaryDofIndices_.assign(numDofs(), false);
273
274 numScv_ = 0;
275 numScvf_ = 0;
276 numBoundaryScvf_ = 0;
277 // Build the SCV and SCV faces
278 for (const auto& element : elements(this->gridView()))
279 {
280 // fill the element map with seeds
281 const auto eIdx = this->elementMapper().index(element);
282
283 // count
284 numScv_ += element.subEntities(dim);
285 numScvf_ += element.subEntities(dim-1);
286
287 // get the element geometry
288 auto elementGeometry = element.geometry();
289 const auto refElement = referenceElement(elementGeometry);
290
291 // instantiate the geometry helper
292 GeometryHelper geometryHelper(elementGeometry);
293
294 // construct the sub control volumes
295 cache_.scvs_[eIdx].resize(elementGeometry.corners());
296 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
297 {
298 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
299
300 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(
301 geometryHelper.getScvCorners(scvLocalIdx),
302 scvLocalIdx,
303 eIdx,
304 dofIdxGlobal
305 );
306 }
307
308 // construct the sub control volume faces
309 LocalIndexType scvfLocalIdx = 0;
310 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
311 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
312 {
313 // find the global and local scv indices this scvf is belonging to
314 std::array<LocalIndexType, 2> localScvIndices{{
315 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
316 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
317 }};
318
319 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
320 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(
321 corners,
322 geometryHelper.normal(corners, localScvIndices),
323 element,
324 scvfLocalIdx,
325 std::move(localScvIndices)
326 );
327 }
328
329 // construct the sub control volume faces on the domain boundary
330 for (const auto& intersection : intersections(this->gridView(), element))
331 {
332 if (intersection.boundary() && !intersection.neighbor())
333 {
334 const auto isGeometry = intersection.geometry();
335 cache_.hasBoundaryScvf_[eIdx] = true;
336
337 // count
338 numScvf_ += isGeometry.corners();
339 numBoundaryScvf_ += isGeometry.corners();
340
341 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
342 {
343 // find the scvs this scvf is belonging to
344 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
345 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
346
347 cache_.scvfs_[eIdx].emplace_back(
348 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
349 intersection.centerUnitOuterNormal(),
350 intersection,
351 isScvfLocalIdx,
352 scvfLocalIdx,
353 std::move(localScvIndices)
354 );
355
356 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
357 static_cast<LocalIndexType>(intersection.indexInInside()),
358 static_cast<LocalIndexType>(isScvfLocalIdx)
359 }});
360
361 // increment local counter
362 scvfLocalIdx++;
363 }
364
365 // add all vertices on the intersection to the set of
366 // boundary vertices
367 const auto fIdx = intersection.indexInInside();
368 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
369 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
370 {
371 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
372 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
373 boundaryDofIndices_[vIdxGlobal] = true;
374 }
375 }
376
377 // inform the grid geometry if we have periodic boundaries
378 else if (periodicGridTraits_.isPeriodic(intersection))
379 {
380 this->setPeriodic();
381
382 // find the mapped periodic vertex of all vertices on periodic boundaries
383 const auto fIdx = intersection.indexInInside();
384 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
385 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
386 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
387 {
388 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
389 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
390 const auto vPos = elementGeometry.corner(vIdx);
391
392 const auto& outside = intersection.outside();
393 const auto outsideGeometry = outside.geometry();
394 for (const auto& isOutside : intersections(this->gridView(), outside))
395 {
396 // only check periodic vertices of the periodic neighbor
397 if (periodicGridTraits_.isPeriodic(isOutside))
398 {
399 const auto fIdxOutside = isOutside.indexInInside();
400 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
401 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
402 {
403 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
404 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
405 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
406 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
407 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
408 }
409 }
410 }
411 }
412 }
413 }
414 }
415
416 // error check: periodic boundaries currently don't work for box in parallel
417 if (this->isPeriodic() && this->gridView().comm().size() > 1)
418 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
419 }
420
421 const FeCache feCache_;
422
423 std::size_t numScv_;
424 std::size_t numScvf_;
425 std::size_t numBoundaryScvf_;
426
427 // vertices on the boundary
428 std::vector<bool> boundaryDofIndices_;
429
430 // a map for periodic boundary vertices
431 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
432
433 Cache cache_;
434
436};
437
445template<class Scalar, class GV, class Traits>
446class BoxFVGridGeometry<Scalar, GV, false, Traits>
447: public BaseGridGeometry<GV, Traits>
448{
451 using GridIndexType = typename IndexTraits<GV>::GridIndex;
452
453 static const int dim = GV::dimension;
454 static const int dimWorld = GV::dimensionworld;
455
456 using Element = typename GV::template Codim<0>::Entity;
457 using CoordScalar = typename GV::ctype;
458
459public:
462 static constexpr DiscretizationMethod discMethod{};
463
467 using LocalView = typename Traits::template LocalView<ThisType, false>;
469 using SubControlVolume = typename Traits::SubControlVolume;
471 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
475 using DofMapper = typename Traits::VertexMapper;
477 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
479 using GridView = GV;
483 using ScvQuadratureRule = typename Traits::ScvQuadratureRule;
485 using ScvfQuadratureRule = typename Traits::ScvfQuadratureRule;
486
488 BoxFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg)
489 : ParentType(std::move(gg))
490 , cache_(*this)
491 , periodicGridTraits_(this->gridView().grid())
492 {
493 update_();
494 }
495
498 : BoxFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView))
499 {}
500
503 const DofMapper& dofMapper() const
504 { return this->vertexMapper(); }
505
507 std::size_t numScv() const
508 { return numScv_; }
509
511 std::size_t numScvf() const
512 { return numScvf_; }
513
516 std::size_t numBoundaryScvf() const
517 { return numBoundaryScvf_; }
518
520 std::size_t numDofs() const
521 { return this->vertexMapper().size(); }
522
523
525 void update(const GridView& gridView)
526 {
527 ParentType::update(gridView);
528 update_();
529 }
530
532 void update(GridView&& gridView)
533 {
534 ParentType::update(std::move(gridView));
535 update_();
536 }
537
539 const FeCache& feCache() const
540 { return feCache_; }
541
543 bool dofOnBoundary(GridIndexType dofIdx) const
544 { return boundaryDofIndices_[dofIdx]; }
545
547 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
548 { return periodicDofMap_.count(dofIdx); }
549
551 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
552 { return periodicDofMap_.at(dofIdx); }
553
555 const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const
556 { return periodicDofMap_; }
557
559 friend inline LocalView localView(const BoxFVGridGeometry& gg)
560 { return { gg.cache_ }; }
561
562private:
563
564 class BoxGridGeometryCache
565 {
566 friend class BoxFVGridGeometry;
567 public:
569 using GeometryHelper = Detail::BoxGeometryHelper_t<GV, Traits>;
570
571 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
572 : gridGeometry_(&gg)
573 {}
574
575 const BoxFVGridGeometry& gridGeometry() const
576 { return *gridGeometry_; }
577
578 private:
579 const BoxFVGridGeometry* gridGeometry_;
580 };
581
582public:
585 using Cache = BoxGridGeometryCache;
586
587private:
588
589 void update_()
590 {
591 boundaryDofIndices_.assign(numDofs(), false);
592
593 // save global data on the grid's scvs and scvfs
594 // TODO do we need those information?
595 numScv_ = 0;
596 numScvf_ = 0;
597 numBoundaryScvf_ = 0;
598 for (const auto& element : elements(this->gridView()))
599 {
600 numScv_ += element.subEntities(dim);
601 numScvf_ += element.subEntities(dim-1);
602
603 const auto elementGeometry = element.geometry();
604 const auto refElement = referenceElement(elementGeometry);
605
606 // store the sub control volume face indices on the domain boundary
607 for (const auto& intersection : intersections(this->gridView(), element))
608 {
609 if (intersection.boundary() && !intersection.neighbor())
610 {
611 const auto isGeometry = intersection.geometry();
612 numScvf_ += isGeometry.corners();
613 numBoundaryScvf_ += isGeometry.corners();
614
615 // add all vertices on the intersection to the set of
616 // boundary vertices
617 const auto fIdx = intersection.indexInInside();
618 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
619 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
620 {
621 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
622 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
623 boundaryDofIndices_[vIdxGlobal] = true;
624 }
625 }
626
627 // inform the grid geometry if we have periodic boundaries
628 else if (periodicGridTraits_.isPeriodic(intersection))
629 {
630 this->setPeriodic();
631
632 // find the mapped periodic vertex of all vertices on periodic boundaries
633 const auto fIdx = intersection.indexInInside();
634 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
635 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
636 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
637 {
638 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
639 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
640 const auto vPos = elementGeometry.corner(vIdx);
641
642 const auto& outside = intersection.outside();
643 const auto outsideGeometry = outside.geometry();
644 for (const auto& isOutside : intersections(this->gridView(), outside))
645 {
646 // only check periodic vertices of the periodic neighbor
647 if (periodicGridTraits_.isPeriodic(isOutside))
648 {
649 const auto fIdxOutside = isOutside.indexInInside();
650 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
651 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
652 {
653 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
654 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
655 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
656 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
657 periodicDofMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
658 }
659 }
660 }
661 }
662 }
663 }
664 }
665
666 // error check: periodic boundaries currently don't work for box in parallel
667 if (this->isPeriodic() && this->gridView().comm().size() > 1)
668 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
669 }
670
671 const FeCache feCache_;
672
673 // Information on the global number of geometries
674 // TODO do we need those information?
675 std::size_t numScv_;
676 std::size_t numScvf_;
677 std::size_t numBoundaryScvf_;
678
679 // vertices on the boundary
680 std::vector<bool> boundaryDofIndices_;
681
682 // a map for periodic boundary vertices
683 std::unordered_map<GridIndexType, GridIndexType> periodicDofMap_;
684
685 Cache cache_;
686
688};
689
690} // end namespace Dumux
691
692#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:44
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:448
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:516
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:559
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:547
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:467
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:471
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:497
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:469
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition: discretization/box/fvgridgeometry.hh:483
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:551
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:539
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:585
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:511
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:475
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:465
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:543
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvgridgeometry.hh:485
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:503
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition: discretization/box/fvgridgeometry.hh:481
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:532
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:473
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:520
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:488
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:507
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:525
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:477
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:555
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:92
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:203
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:258
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:169
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:195
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:164
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:151
BasicGridGeometry_t< GV, Traits > BasicGridGeometry
export basic grid geometry type for the alternative constructor
Definition: discretization/box/fvgridgeometry.hh:109
typename Traits::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition: discretization/box/fvgridgeometry.hh:127
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:183
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:160
typename Traits::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvgridgeometry.hh:129
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/box/fvgridgeometry.hh:121
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:115
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:147
const std::unordered_map< GridIndexType, GridIndexType > & periodicDofMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:199
BoxFVGridGeometry(const GridView &gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:141
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:132
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:117
typename PeriodicGridTraits< typename GV::Grid >::SupportsPeriodicity SupportsPeriodicity
export whether the grid(geometry) supports periodicity
Definition: discretization/box/fvgridgeometry.hh:125
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/box/fvgridgeometry.hh:111
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:187
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:155
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:191
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:176
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:113
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/box/fvgridgeometry.hh:119
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:81
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:98
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:63
Quadrature rule traits for discretization schemes.
Definition: quadraturerules.hh:84
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
Definition: periodicgridtraits.hh:24