3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/box/fvgridgeometry.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
27#define DUMUX_DISCRETIZATION_BOX_GRID_FVGEOMETRY_HH
28
29#include <utility>
30#include <unordered_map>
31#include <array>
32#include <vector>
33
34#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
35
45
46namespace Dumux {
47
54template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
56: public MapperTraits
57{
60
61 template<class GridGeometry, bool enableCache>
63};
64
71template<class Scalar,
72 class GridView,
73 bool enableGridGeometryCache = false,
76
83template<class Scalar, class GV, class Traits>
84class BoxFVGridGeometry<Scalar, GV, true, Traits>
85: public BaseGridGeometry<GV, Traits>
86{
89 using GridIndexType = typename IndexTraits<GV>::GridIndex;
90 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
91
92 using Element = typename GV::template Codim<0>::Entity;
93 using CoordScalar = typename GV::ctype;
94 static const int dim = GV::dimension;
95 static const int dimWorld = GV::dimensionworld;
96
97 using GeometryHelper = BoxGeometryHelper<GV, dim,
98 typename Traits::SubControlVolume,
99 typename Traits::SubControlVolumeFace>;
100
101 // Remove this after release 3.6 when the deprecated local view constructor is removed
102 using LV_ = typename Traits::template LocalView<ThisType, true>;
103 friend LV_;
104
105public:
108 static constexpr DiscretizationMethod discMethod{};
109
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;
124
127 : ParentType(gridView)
128 , cache_(*this)
129 {
130 update_();
131 }
132
135 const DofMapper& dofMapper() const
136 { return this->vertexMapper(); }
137
139 std::size_t numScv() const
140 { return numScv_; }
141
143 std::size_t numScvf() const
144 { return numScvf_; }
145
148 std::size_t numBoundaryScvf() const
149 { return numBoundaryScvf_; }
150
152 std::size_t numDofs() const
153 { return this->vertexMapper().size(); }
154
155
157 void update(const GridView& gridView)
158 {
159 ParentType::update(gridView);
160 update_();
161 }
162
164 void update(GridView&& gridView)
165 {
166 ParentType::update(std::move(gridView));
167 update_();
168 }
169
171 const FeCache& feCache() const
172 { return feCache_; }
173
175 [[deprecated("Will be removed after release 3.6")]]
176 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const
177 { return cache_.scvs_[eIdx]; }
178
180 [[deprecated("Will be removed after release 3.6")]]
181 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const
182 { return cache_.scvfs_[eIdx]; }
183
185 bool dofOnBoundary(GridIndexType dofIdx) const
186 { return boundaryDofIndices_[dofIdx]; }
187
189 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
190 { return periodicVertexMap_.count(dofIdx); }
191
193 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
194 { return periodicVertexMap_.at(dofIdx); }
195
197 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
198 { return periodicVertexMap_; }
199
201 [[deprecated("Will be removed after release 3.6")]]
202 bool hasBoundaryScvf(GridIndexType eIdx) const
203 { return cache_.hasBoundaryScvf_[eIdx]; }
204
206 friend inline LocalView localView(const BoxFVGridGeometry& gg)
207 { return { gg.cache_ }; }
208
209private:
210
211 class BoxGridGeometryCache
212 {
213 friend class BoxFVGridGeometry;
214 public:
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 void update_()
262 {
263 cache_.clear_();
264
265 const auto numElements = this->gridView().size(0);
266 cache_.scvs_.resize(numElements);
267 cache_.scvfs_.resize(numElements);
268 cache_.hasBoundaryScvf_.resize(numElements, false);
269
270 boundaryDofIndices_.assign(numDofs(), false);
271
272 numScv_ = 0;
273 numScvf_ = 0;
274 numBoundaryScvf_ = 0;
275 // Build the SCV and SCV faces
276 for (const auto& element : elements(this->gridView()))
277 {
278 // fill the element map with seeds
279 const auto eIdx = this->elementMapper().index(element);
280
281 // count
282 numScv_ += element.subEntities(dim);
283 numScvf_ += element.subEntities(dim-1);
284
285 // get the element geometry
286 auto elementGeometry = element.geometry();
287 const auto refElement = referenceElement(elementGeometry);
288
289 // instantiate the geometry helper
290 GeometryHelper geometryHelper(elementGeometry);
291
292 // construct the sub control volumes
293 cache_.scvs_[eIdx].resize(elementGeometry.corners());
294 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
295 {
296 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
297
298 cache_.scvs_[eIdx][scvLocalIdx] = SubControlVolume(
299 geometryHelper,
300 scvLocalIdx,
301 eIdx,
302 dofIdxGlobal
303 );
304 }
305
306 // construct the sub control volume faces
307 LocalIndexType scvfLocalIdx = 0;
308 cache_.scvfs_[eIdx].resize(element.subEntities(dim-1));
309 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
310 {
311 // find the global and local scv indices this scvf is belonging to
312 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
313 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
314
315 cache_.scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(
316 geometryHelper,
317 element,
318 elementGeometry,
319 scvfLocalIdx,
320 std::move(localScvIndices),
321 false
322 );
323 }
324
325 // construct the sub control volume faces on the domain boundary
326 for (const auto& intersection : intersections(this->gridView(), element))
327 {
328 if (intersection.boundary() && !intersection.neighbor())
329 {
330 const auto isGeometry = intersection.geometry();
331 cache_.hasBoundaryScvf_[eIdx] = true;
332
333 // count
334 numScvf_ += isGeometry.corners();
335 numBoundaryScvf_ += isGeometry.corners();
336
337 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
338 {
339 // find the scvs this scvf is belonging to
340 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
341 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
342
343 cache_.scvfs_[eIdx].emplace_back(
344 geometryHelper,
345 intersection,
346 isGeometry,
347 isScvfLocalIdx,
348 scvfLocalIdx,
349 std::move(localScvIndices),
350 true
351 );
352
353 cache_.scvfBoundaryGeometryKeys_[eIdx].emplace_back(std::array<LocalIndexType, 2>{{
354 static_cast<LocalIndexType>(intersection.indexInInside()),
355 static_cast<LocalIndexType>(isScvfLocalIdx)
356 }});
357
358 // increment local counter
359 scvfLocalIdx++;
360 }
361
362 // add all vertices on the intersection to the set of
363 // boundary vertices
364 const auto fIdx = intersection.indexInInside();
365 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
366 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
367 {
368 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
369 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
370 boundaryDofIndices_[vIdxGlobal] = true;
371 }
372 }
373
374 // inform the grid geometry if we have periodic boundaries
375 else if (intersection.boundary() && intersection.neighbor())
376 {
377 this->setPeriodic();
378
379 // find the mapped periodic vertex of all vertices on periodic boundaries
380 const auto fIdx = intersection.indexInInside();
381 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
382 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
383 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
384 {
385 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
386 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
387 const auto vPos = elementGeometry.corner(vIdx);
388
389 const auto& outside = intersection.outside();
390 const auto outsideGeometry = outside.geometry();
391 for (const auto& isOutside : intersections(this->gridView(), outside))
392 {
393 // only check periodic vertices of the periodic neighbor
394 if (isOutside.boundary() && isOutside.neighbor())
395 {
396 const auto fIdxOutside = isOutside.indexInInside();
397 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
398 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
399 {
400 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
401 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
402 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
403 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
404 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
405 }
406 }
407 }
408 }
409 }
410 }
411 }
412
413 // error check: periodic boundaries currently don't work for box in parallel
414 if (this->isPeriodic() && this->gridView().comm().size() > 1)
415 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
416 }
417
418 const FeCache feCache_;
419
420 std::size_t numScv_;
421 std::size_t numScvf_;
422 std::size_t numBoundaryScvf_;
423
424 // vertices on the boundary
425 std::vector<bool> boundaryDofIndices_;
426
427 // a map for periodic boundary vertices
428 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
429
430 Cache cache_;
431};
432
440template<class Scalar, class GV, class Traits>
441class BoxFVGridGeometry<Scalar, GV, false, Traits>
442: public BaseGridGeometry<GV, Traits>
443{
446 using GridIndexType = typename IndexTraits<GV>::GridIndex;
447
448 static const int dim = GV::dimension;
449 static const int dimWorld = GV::dimensionworld;
450
451 using Element = typename GV::template Codim<0>::Entity;
452 using CoordScalar = typename GV::ctype;
453
454 // Remove this after release 3.6 when the deprecated local view constructor is removed
455 using LV_ = typename Traits::template LocalView<ThisType, false>;
456 friend LV_;
457
458public:
461 static constexpr DiscretizationMethod discMethod{};
462
464 using LocalView = typename Traits::template LocalView<ThisType, false>;
466 using SubControlVolume = typename Traits::SubControlVolume;
468 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
472 using DofMapper = typename Traits::VertexMapper;
474 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
476 using GridView = GV;
477
480 : ParentType(gridView)
481 , cache_(*this)
482 {
483 update_();
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 periodicVertexMap_.count(dofIdx); }
534
536 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
537 { return periodicVertexMap_.at(dofIdx); }
538
540 const std::unordered_map<GridIndexType, GridIndexType>& periodicVertexMap() const
541 { return periodicVertexMap_; }
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:
553 explicit BoxGridGeometryCache(const BoxFVGridGeometry& gg)
554 : gridGeometry_(&gg)
555 {}
556
557 const BoxFVGridGeometry& gridGeometry() const
558 { return *gridGeometry_; }
559
560 private:
561 const BoxFVGridGeometry* gridGeometry_;
562 };
563
564public:
567 using Cache = BoxGridGeometryCache;
568
569private:
570
571 void update_()
572 {
573 boundaryDofIndices_.assign(numDofs(), false);
574
575 // save global data on the grid's scvs and scvfs
576 // TODO do we need those information?
577 numScv_ = 0;
578 numScvf_ = 0;
579 numBoundaryScvf_ = 0;
580 for (const auto& element : elements(this->gridView()))
581 {
582 numScv_ += element.subEntities(dim);
583 numScvf_ += element.subEntities(dim-1);
584
585 const auto elementGeometry = element.geometry();
586 const auto refElement = referenceElement(elementGeometry);
587
588 // store the sub control volume face indices on the domain boundary
589 for (const auto& intersection : intersections(this->gridView(), element))
590 {
591 if (intersection.boundary() && !intersection.neighbor())
592 {
593 const auto isGeometry = intersection.geometry();
594 numScvf_ += isGeometry.corners();
595 numBoundaryScvf_ += isGeometry.corners();
596
597 // add all vertices on the intersection to the set of
598 // boundary vertices
599 const auto fIdx = intersection.indexInInside();
600 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
601 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
602 {
603 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
604 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
605 boundaryDofIndices_[vIdxGlobal] = true;
606 }
607 }
608
609 // inform the grid geometry if we have periodic boundaries
610 else if (intersection.boundary() && intersection.neighbor())
611 {
612 this->setPeriodic();
613
614 // find the mapped periodic vertex of all vertices on periodic boundaries
615 const auto fIdx = intersection.indexInInside();
616 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
617 const auto eps = 1e-7*(elementGeometry.corner(1) - elementGeometry.corner(0)).two_norm();
618 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
619 {
620 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
621 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
622 const auto vPos = elementGeometry.corner(vIdx);
623
624 const auto& outside = intersection.outside();
625 const auto outsideGeometry = outside.geometry();
626 for (const auto& isOutside : intersections(this->gridView(), outside))
627 {
628 // only check periodic vertices of the periodic neighbor
629 if (isOutside.boundary() && isOutside.neighbor())
630 {
631 const auto fIdxOutside = isOutside.indexInInside();
632 const auto numFaceVertsOutside = refElement.size(fIdxOutside, 1, dim);
633 for (int localVIdxOutside = 0; localVIdxOutside < numFaceVertsOutside; ++localVIdxOutside)
634 {
635 const auto vIdxOutside = refElement.subEntity(fIdxOutside, 1, localVIdxOutside, dim);
636 const auto vPosOutside = outsideGeometry.corner(vIdxOutside);
637 const auto shift = std::abs((this->bBoxMax()-this->bBoxMin())*intersection.centerUnitOuterNormal());
638 if (std::abs((vPosOutside-vPos).two_norm() - shift) < eps)
639 periodicVertexMap_[vIdxGlobal] = this->vertexMapper().subIndex(outside, vIdxOutside, dim);
640 }
641 }
642 }
643 }
644 }
645 }
646 }
647
648 // error check: periodic boundaries currently don't work for box in parallel
649 if (this->isPeriodic() && this->gridView().comm().size() > 1)
650 DUNE_THROW(Dune::NotImplemented, "Periodic boundaries for box method for parallel simulations!");
651 }
652
653 const FeCache feCache_;
654
655 // Information on the global number of geometries
656 // TODO do we need those information?
657 std::size_t numScv_;
658 std::size_t numScvf_;
659 std::size_t numBoundaryScvf_;
660
661 // vertices on the boundary
662 std::vector<bool> boundaryDofIndices_;
663
664 // a map for periodic boundary vertices
665 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap_;
666
667 Cache cache_;
668};
669
670} // end namespace Dumux
671
672#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for grid geometries.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all grid geometries.
Definition: basegridgeometry.hh:61
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:69
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:53
The default traits for the box finite volume grid geometry Defines the scv and scvf types and the map...
Definition: discretization/box/fvgridgeometry.hh:57
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:75
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: discretization/box/fvgridgeometry.hh:86
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: discretization/box/fvgridgeometry.hh:176
friend LocalView localView(const BoxFVGridGeometry &gg)
local view of this object (constructed with the internal cache)
Definition: discretization/box/fvgridgeometry.hh:206
BoxGridGeometryCache Cache
Definition: discretization/box/fvgridgeometry.hh:258
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: discretization/box/fvgridgeometry.hh:181
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:197
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:157
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:193
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvgridgeometry.hh:202
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:152
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvgridgeometry.hh:139
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/box/fvgridgeometry.hh:171
std::size_t numBoundaryScvf() const
Definition: discretization/box/fvgridgeometry.hh:148
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:135
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/box/fvgridgeometry.hh:117
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
BoxFVGridGeometry(const GridView gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:126
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:185
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/box/fvgridgeometry.hh:143
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary.
Definition: discretization/box/fvgridgeometry.hh:189
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/box/fvgridgeometry.hh:164
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh: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:443
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:464
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:468
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/box/fvgridgeometry.hh:466
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:567
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:472
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/box/fvgridgeometry.hh:528
BoxFVGridGeometry(const GridView gridView)
Constructor.
Definition: discretization/box/fvgridgeometry.hh:479
const DofMapper & dofMapper() const
Definition: discretization/box/fvgridgeometry.hh:488
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:470
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/box/fvgridgeometry.hh:505
const std::unordered_map< GridIndexType, GridIndexType > & periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/box/fvgridgeometry.hh:540
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:474
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:72
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:74
Definition: method.hh:58
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.